GRAVI Pipeline Reference Manual 1.9.6
Loading...
Searching...
No Matches
gravi_dfs.c
Go to the documentation of this file.
1/* $Id: gravi_dfs.c,v 1.6 2011/04/31 06:10:40 llundin Exp $
2 *
3 * This file is part of the GRAVI Pipeline
4 * Copyright (C) 2002,2003 European Southern Observatory
5 *
6 * This program is free software; you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License as published by
8 * the Free Software Foundation; either version 2 of the License, or
9 * (at your option) any later version.
10 *
11 * This program 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 General Public License for more details.
15 *
16 * You should have received a copy of the GNU General Public License
17 * along with this program; if not, write to the Free Software
18 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19 */
20
35/*
36 * History :
37 * 12/11/2018 add gravi_frameset_extract_static_param
38 * 04/12/2018 add gravi_frameset_extract_wave_param
39 */
40
41/*-----------------------------------------------------------------------------
42 Includes
43 -----------------------------------------------------------------------------*/
44
45#ifdef HAVE_CONFIG_H
46#include <config.h>
47#endif
48
49#include "time.h"
50#include <string.h>
51#include <math.h>
52#include <cpl.h>
53
54#include "gravi_dfs.h"
55#include "gravi_utils.h"
56
57/*-----------------------------------------------------------------------------
58 Functions code
59 -----------------------------------------------------------------------------*/
60
62{
63 cpl_msg_info(__func__, "**********************************************");
64 cpl_msg_info(__func__, " Welcome to GRAVITY Pipeline release %s",
65 PACKAGE_VERSION);
66 cpl_msg_info(__func__, " Last rebuilt at %s %s",__DATE__,__TIME__);
67 cpl_msg_info(__func__, "**********************************************");
68}
69
70/*----------------------------------------------------------------------------*/
76/*----------------------------------------------------------------------------*/
77
78cpl_error_code gravi_dfs_set_groups(cpl_frameset *set) {
79 cpl_ensure_code(set, CPL_ERROR_NULL_INPUT);
80
81 cpl_errorstate prestate = cpl_errorstate_get();
82 cpl_frame *frame = NULL;
83 int i, nb_frame;
84 nb_frame = cpl_frameset_get_size(set);
85
86 /* Loop on frames */
87 for (i = 0; i < nb_frame; i++) {
88
89 frame = cpl_frameset_get_position(set, i);
90 const char *tag = cpl_frame_get_tag(frame);
91
92 if (tag == NULL) {
93 cpl_msg_warning(cpl_func, "Frame %d has no tag", i);
94 } else if ((!strcmp(tag, GRAVI_DARK_RAW)) ||
95 (!strcmp(tag, GRAVI_FLAT_RAW)) ||
96 (!strcmp(tag, GRAVI_WAVE_RAW)) ||
97 (!strcmp(tag, GRAVI_WAVESC_RAW)) ||
98 (!strcmp(tag, GRAVI_WAVELAMP_RAW)) ||
99 (!strcmp(tag, GRAVI_P2VM_RAW)) ||
100 (!strcmp(tag, GRAVI_PIEZOTF_RAW)) ||
101 (!strcmp(tag, GRAVI_SINGLE_SCIENCE_RAW)) ||
102 (!strcmp(tag, GRAVI_SINGLE_CALIB_RAW)) ||
103 (!strcmp(tag, GRAVI_DUAL_SCIENCE_RAW)) ||
104 (!strcmp(tag, GRAVI_DUAL_CALIB_RAW)) ||
105 (!strcmp(tag, GRAVI_DUAL_SKY_RAW)) ||
106 (!strcmp(tag, GRAVI_SINGLE_SKY_RAW)) ||
107 (!strcmp(tag, GRAVI_VIS_SINGLE_CALIB)) ||
108 (!strcmp(tag, GRAVI_VISPHI_SINGLE_CALIB)) ||
109 (!strcmp(tag, GRAVI_VIS_SINGLE_SCIENCE)) ||
110 (!strcmp(tag, GRAVI_VIS_DUAL_CALIB)) ||
111 (!strcmp(tag, GRAVI_VISPHI_DUAL_CALIB)) ||
112 (!strcmp(tag, GRAVI_VIS_DUAL_SCIENCE)) ||
113 (!strcmp(tag, GRAVI_P2VMRED_SINGLE_CALIB)) ||
114 (!strcmp(tag, GRAVI_P2VMRED_SINGLE_SCIENCE)) ||
115 (!strcmp(tag, GRAVI_P2VMRED_DUAL_CALIB)) ||
116 (!strcmp(tag, GRAVI_P2VMRED_DUAL_SCIENCE)) ||
117 (!strcmp(tag, GRAVI_MIRA_INPUT_PROCATG)) ||
118 (!strcmp(tag, GRAVI_VIS_SINGLE_CALIBRATED)) ||
119 (!strcmp(tag, GRAVI_VIS_DUAL_CALIBRATED)) ||
120 (!strcmp(tag, GRAVI_DISP_RAW))) {
121 /* RAW frames */
122 cpl_frame_set_group(frame, CPL_FRAME_GROUP_RAW);
123 } else if ((!strcmp(tag, GRAVI_DARK_MAP)) ||
124 (!strcmp(tag, GRAVI_FLAT_MAP)) ||
125 (!strcmp(tag, GRAVI_WAVE_MAP)) ||
126 (!strcmp(tag, GRAVI_P2VM_MAP)) ||
127 (!strcmp(tag, GRAVI_BAD_MAP)) ||
128 (!strcmp(tag, GRAVI_BIASMASK_MAP)) ||
129 (!strcmp(tag, GRAVI_PIEZOTF_MAP)) ||
130 (!strcmp(tag, GRAVI_PREPROC)) ||
131 (!strcmp(tag, GRAVI_TF_SINGLE_SCIENCE)) ||
132 (!strcmp(tag, GRAVI_TF_SINGLE_CALIB)) ||
133 (!strcmp(tag, GRAVI_VISPHI_TF_SINGLE_CALIB)) ||
134 (!strcmp(tag, GRAVI_WAVELAMP_MAP)) ||
135 (!strcmp(tag, GRAVI_TF_DUAL_SCIENCE)) ||
136 (!strcmp(tag, GRAVI_TF_DUAL_CALIB)) ||
137 (!strcmp(tag, GRAVI_VISPHI_TF_DUAL_CALIB)) ||
138 (!strcmp(tag, GRAVI_ZP_CAL)) ||
139 (!strcmp(tag, GRAVI_DISP_VIS)) ||
140 (!strcmp(tag, GRAVI_DIAMETER_CAT)) ||
141 (!strcmp(tag, GRAVI_DISP_MODEL)) ||
142 (!strcmp(tag, GRAVI_DIODE_POSITION)) ||
143 (!strcmp(tag, GRAVI_KEY_PATCH)) ||
144 (!strcmp(tag, GRAVI_ASTRO_CAL_PHASEREF)) ||
145 (!strcmp(tag, GRAVI_ASTRO_TARGET)) ||
146 (!strcmp(tag, GRAVI_ASTRO_SWAP))) {
147 /* CALIB frames */
148 cpl_frame_set_group(frame, CPL_FRAME_GROUP_CALIB);
149 } else if ((!strcmp(tag, GRAVI_MIRA_OUTPUT_PROCATG)) ||
150 (!strcmp(tag, GRAVI_NAB_CAL))) {
151 /* PRODUCT frames */
152 cpl_frame_set_group(frame, CPL_FRAME_GROUP_PRODUCT);
153 }
154 }
155
156 if (!cpl_errorstate_is_equal(prestate)) {
157 return cpl_error_set_message(cpl_func, cpl_error_get_code(),
158 "Could not identify RAW and CALIB "
159 "frames");
160 }
161
162 return CPL_ERROR_NONE;
163}
164/*----------------------------------------------------------------------------*/
171/*----------------------------------------------------------------------------*/
172
173cpl_error_code gravi_parameter_disable (cpl_parameter * p)
174{
175 cpl_ensure_code (p, CPL_ERROR_NULL_INPUT);
176 cpl_parameter_disable (p, CPL_PARAMETER_MODE_CLI);
177 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
178 cpl_parameter_disable (p, CPL_PARAMETER_MODE_CFG);
179 return CPL_ERROR_NONE;
180}
181
182/*----------------------------------------------------------------------------*/
189/*----------------------------------------------------------------------------*/
190
191cpl_parameter * gravi_parameter_add_badpix (cpl_parameterlist *self)
192{
193 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
194
195 cpl_parameter *p;
196 p = cpl_parameter_new_value ("gravity.calib.bad-dark-threshold", CPL_TYPE_INT,
197 "the rms factor for "
198 "dark bad pixel threshold",
199 "gravity.calib", 10);
200 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "bad-dark-threshold");
201 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
202 cpl_parameterlist_append (self, p);
203
204 p = cpl_parameter_new_value ("gravity.calib.lowflux-pixels-ft", CPL_TYPE_BOOL,
205 "Flag as bad pixels all pixels on the FT which "
206 "have a non-sgnificant flux "
207 "(less than 33% of neighbouring pixel)"
208 "to increase SNR on faint targets",
209 "gravity.calib", FALSE);
210 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "lowflux-pixels-ft");
211 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
212 cpl_parameterlist_append (self, p);
213
214 p = cpl_parameter_new_value ("gravity.calib.bad-pixel-A-ft", CPL_TYPE_INT,
215 "flag a given pixel as bad on the FT detector"
216 "value is position of pixel (starting at 1)"
217 "a position of 0 means no pixel is flagged",
218 "gravity.calib", 0);
219 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "bad-pixel-A-ft");
220 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
221 cpl_parameterlist_append (self, p);
222
223 p = cpl_parameter_new_value ("gravity.calib.bad-pixel-B-ft", CPL_TYPE_INT,
224 "flag a second pixel as bad on the FT detector"
225 "value is position of pixel (starting at 1)"
226 "a position of 0 means no pixel is flagged",
227 "gravity.calib", 0);
228 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "bad-pixel-B-ft");
229 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
230 cpl_parameterlist_append (self, p);
231
232 return p;
233}
234
235/*----------------------------------------------------------------------------*/
242/*----------------------------------------------------------------------------*/
243
244cpl_parameter * gravi_parameter_add_pcacalib (cpl_parameterlist *self)
245{
246 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
247
248 cpl_parameter *p;
249
250 /* Number of PCA components */
251 p = cpl_parameter_new_value("gravity.calib.pca-components", CPL_TYPE_INT,
252 "The number of PCA components to compute.",
253 "gravity.calib", 1);
254 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "pca-components");
255 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
256 cpl_parameterlist_append (self, p);
257
258 /* Tracking ratio threshold */
259 p = cpl_parameter_new_value("gravity.calib.pca-tracking-ratio", CPL_TYPE_INT,
260 "The minimum tracking ratio to accept calibration frames for.",
261 "gravity.calib", 90);
262 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "pca-tracking-ratio");
263 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
264 cpl_parameterlist_append (self, p);
265
266 /* Outlier cleaning window size */
267 p = cpl_parameter_new_value("gravity.calib.pca-clean-size", CPL_TYPE_INT,
268 "The window size to use for outlier cleaning.",
269 "gravity.calib", 20);
270 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "pca-clean-size");
271 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
272 cpl_parameterlist_append (self, p);
273
274 /* Tracking ratio threshold */
275 p = cpl_parameter_new_value("gravity.calib.pca-clean-nstd", CPL_TYPE_DOUBLE,
276 "The sigma-clip n_std for outlier cleaning.",
277 "gravity.calib", 5.0);
278 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "pca-clean-nstd");
279 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
280 cpl_parameterlist_append (self, p);
281
282 /* method for fitting PCA components */
283 const char *PCA_FIT_TYPES[2] = {"POLYNOMIAL", "SPLINE"};
284 p = cpl_parameter_new_enum("gravity.calib.pca-fit-type", CPL_TYPE_STRING,
285 "The method to use for fitting the PCA components.",
286 "gravity.calib",
287 PCA_FIT_TYPES[1], 2, PCA_FIT_TYPES[0], PCA_FIT_TYPES[1]
288 );
289 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "pca-fit-type");
290 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
291 cpl_parameterlist_append (self, p);
292
293 /* degree of model for fitting PCA components */
294 p = cpl_parameter_new_value("gravity.calib.pca-fit-degree", CPL_TYPE_INT,
295 "The polynomial fit degree, or number of spline coefficients.",
296 "gravity.calib", 20);
297 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "pca-fit-degree");
298 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
299 cpl_parameterlist_append (self, p);
300
301 /* whether to output the phase residuals for inspection */
302 p = cpl_parameter_new_value("gravity.calib.pca-save-residuals", CPL_TYPE_BOOL,
303 "Also save the residuals from the PCA fitting for inspection.",
304 "gravity.calib", CPL_FALSE);
305 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "pca-save-residuals");
306 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
307 cpl_parameterlist_append (self, p);
308
309 return p;
310}
311
312/*----------------------------------------------------------------------------*/
319/*----------------------------------------------------------------------------*/
320
321cpl_parameter * gravi_parameter_add_pca (cpl_parameterlist *self)
322{
323 /* Whether to use PCA method to remove VISPHI systematics */
324 cpl_parameter *p = cpl_parameter_new_value ("gravity.vis.flatten-visphi", CPL_TYPE_BOOL,
325 "Use the PCA calibrator to flatten the VISPHI.",
326 "gravity.vis", CPL_FALSE);
327 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "flatten-visphi");
328 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
329 cpl_parameterlist_append (self, p);
330
331 return p;
332}
333
334/*----------------------------------------------------------------------------*/
341/*----------------------------------------------------------------------------*/
342
343cpl_parameter * gravi_parameter_add_profile (cpl_parameterlist *self)
344{
345 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
346
347 cpl_parameter *p;
348
349 /* Method for profile */
350 p = cpl_parameter_new_enum ("gravity.calib.profile-mode", CPL_TYPE_STRING,
351 "Method to compute the extraction profile. "
352 "PROFILE corresponds to the pixel intensities measured in the "
353 "FLAT files (Gaussian like with FWHM of approx 1.5 pixel). "
354 "This is the AUTO option for the Low and Med spectral resolution. "
355 "GAUSS corresponds to a Gaussian fit of the (non-zero) pixel intensities measured "
356 "in the FLAT files. BOX corresponds to a box-card of 6 pixels centered "
357 "on the spectra measured in the FLAT files. This is the AUTO option for High "
358 "spectral resolution",
359 "gravity.calib", "AUTO",
360 4, "AUTO", "PROFILE", "GAUSS", "BOX");
361 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "profile-mode");
362 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
363 cpl_parameterlist_append (self, p);
364
365 /* How to deal with bad-pixels */
366 p = cpl_parameter_new_value ("gravity.calib.force-badpix-to-zero", CPL_TYPE_BOOL,
367 "Force the badpixel to zero in profile",
368 "gravity.calib", TRUE);
369 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "force-badpix-to-zero");
370 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
371 cpl_parameterlist_append (self, p);
372
373 /* The width of the profile element */
374 p = cpl_parameter_new_value ("gravity.calib.profile-width", CPL_TYPE_INT,
375 "Width of the detector window extracted around the default "
376 "position of each spectrum, and on which the profile "
377 "will be applied to perform the extraction.",
378 "gravity.calib", 6);
379 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "profile-width");
380 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
381 cpl_parameterlist_append (self, p);
382
383 return p;
384}
385
386/*----------------------------------------------------------------------------*/
393/*----------------------------------------------------------------------------*/
394cpl_parameter * gravi_parameter_add_preproc (cpl_parameterlist *self)
395{
396 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
397
398 cpl_parameter *p = NULL;
399
400 /* Method for interpolation */
401// p = cpl_parameter_new_value ("gravity.preproc.interp-3pixels", CPL_TYPE_BOOL,
402// "Interpolate with 3 pixels",
403// "gravity.calib", FALSE);
404// cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "interp-3pixels");
405// cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
406// cpl_parameterlist_append (self, p);
407
408 return p;
409}
410
411
412/*----------------------------------------------------------------------------*/
419/*----------------------------------------------------------------------------*/
420cpl_parameter * gravi_parameter_add_wave (cpl_parameterlist *self)
421{
422 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
423
424 cpl_parameter *p;
425
426 /* Method for profile */
427 p = cpl_parameter_new_value ("gravity.calib.force-wave-ft-equal", CPL_TYPE_BOOL,
428 "Force the spatial order of the wavelength 2D fit for FT to "
429 "zero (so all region share the same calibration). "
430 "This is used to build the P2VM calibration of the TAC "
431 "real-time code running on the instrument ifself.",
432 "gravity.calib", FALSE);
433 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "force-wave-ft-equal");
434 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
435 cpl_parameterlist_append (self, p);
436
437 /* Spectral order of 2D wave fit */
438// p = cpl_parameter_new_value ("gravity.calib.wave-spectral-order", CPL_TYPE_INT,
439// "Set the spatial order of the wavelength 2D fit for SC.",
440// "gravity.calib", 3);
441// cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "wave-spectral-order");
442// cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
443// cpl_parameterlist_append (self, p);
444
445 /* Wave determination */
446// p = cpl_parameter_new_enum ("gravity.calib.wave-mode", CPL_TYPE_STRING,
447// "Define the way the wavelength are computed.\n "
448// "PIXEL to compute the wavelength per pixels\n "
449// "BASELINE to compute the wavelength per baseline (ABCD)\n "
450// "AUTO to compute the wavelength per pixels in LOW, and per baseline otherwise",
451// "gravity.calib", "AUTO", 3,
452// "PIXEL","BASELINE","AUTO");
453// cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "wave-mode");
454// cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
455// cpl_parameterlist_append (self, p);
456
457 return p;
458}
459
460
461cpl_parameter * gravi_parameter_add_static_name (cpl_parameterlist *self)
462{
463 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
464
465 cpl_parameter *p;
466 p = cpl_parameter_new_value ("gravity.dfs.static-name", CPL_TYPE_BOOL,
467 "Use static names for the products (for ESO)",
468 "gravity.dfs", FALSE);
469 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "static-name");
470 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
471 cpl_parameterlist_append (self, p);
472
473 return p;
474}
475
476cpl_parameter * gravi_parameter_add_debug_file (cpl_parameterlist *self)
477{
478 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
479
480 cpl_parameter *p;
481 p = cpl_parameter_new_value ("gravity.dfs.debug-file", CPL_TYPE_BOOL,
482 "Save additional debug file(s)",
483 "gravity.dfs", FALSE);
484 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "debug-file");
485 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
486 cpl_parameterlist_append (self, p);
487
488 return p;
489}
490
491cpl_parameter * gravi_parameter_add_biassub_file (cpl_parameterlist *self)
492{
493 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
494
495 cpl_parameter *p;
496 p = cpl_parameter_new_value ("gravity.dfs.bias-subtracted-file", CPL_TYPE_BOOL,
497 "Save the BIAS_SUBTRACTED intermediate product",
498 "gravity.dfs", FALSE);
499 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "bias-subtracted-file");
500 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
501 cpl_parameterlist_append (self, p);
502
503 return p;
504}
505
506cpl_parameter * gravi_parameter_add_spectrum_file (cpl_parameterlist *self)
507{
508 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
509
510 cpl_parameter *p;
511 p = cpl_parameter_new_value ("gravity.dfs.spectrum-file", CPL_TYPE_BOOL,
512 "Save the SPECTRUM intermediate product",
513 "gravity.dfs", FALSE);
514 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "spectrum-file");
515 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
516 cpl_parameterlist_append (self, p);
517
518 return p;
519}
520
521cpl_parameter * gravi_parameter_add_preproc_file (cpl_parameterlist *self)
522{
523 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
524
525 cpl_parameter *p;
526 p = cpl_parameter_new_value ("gravity.dfs.preproc-file", CPL_TYPE_BOOL,
527 "Save the PREPROC intermediate product",
528 "gravity.dfs", FALSE);
529 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "preproc-file");
530 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
531 cpl_parameterlist_append (self, p);
532
533 return p;
534}
535
536cpl_parameter * gravi_parameter_add_p2vmred_file (cpl_parameterlist *self)
537{
538 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
539
540 cpl_parameter *p;
541 p = cpl_parameter_new_value ("gravity.dfs.p2vmred-file", CPL_TYPE_BOOL,
542 "Save the P2VMRED intermediate product",
543 "gravity.dfs", FALSE);
544 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "p2vmreduced-file");
545 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
546 cpl_parameterlist_append (self, p);
547
548 return p;
549}
550
551cpl_parameter * gravi_parameter_add_vis_file (cpl_parameterlist *self)
552{
553 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
554
555 cpl_parameter *p;
556 p = cpl_parameter_new_value ("gravity.dfs.vis-file", CPL_TYPE_BOOL,
557 "Save the VIS intermediate product",
558 "gravity.dfs", FALSE);
559 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "vis-file");
560 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
561 cpl_parameterlist_append (self, p);
562
563 return p;
564}
565
566cpl_parameter * gravi_parameter_add_astro_file (cpl_parameterlist *self)
567{
568 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
569
570 cpl_parameter *p;
571 p = cpl_parameter_new_value ("gravity.dfs.astro-file", CPL_TYPE_BOOL,
572 "Save the ASTROREDUCED intermediate product",
573 "gravity.dfs", FALSE);
574 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "astro-file");
575 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
576 cpl_parameterlist_append (self, p);
577
578 return p;
579}
580
581cpl_parameter * gravi_parameter_add_biasmethod (cpl_parameterlist *self)
582{
583 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
584
585 cpl_parameter *p;
586 p = cpl_parameter_new_enum ("gravity.preproc.bias-method", CPL_TYPE_STRING,
587 "Method to average the biaspixels when cleaning-up\n "
588 "the SC detector (only applied to MED and LOW). Ideally\n "
589 "the same value shall be used when reducing the DARK\n "
590 "with gravity_dark and the OBJECT with gravity_vis. \n"
591 "AUTO is equivalent to MASKED_MEDIAN_PER_COLUMN if the data\n"
592 "contains in the IMAGING_DETECTOR_SC extension the\n"
593 "LEFT, HALFLEFT, CENTER, HALFRIGHT and RIGHT columns.\n"
594 "Otherwise it is like MEDIAN.\n",
595 "gravity.preproc", "AUTO",
596 4, "AUTO", "MEDIAN", "MEDIAN_PER_COLUMN",
597 "MASKED_MEDIAN_PER_COLUMN");
598 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "bias-method");
599 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
600 cpl_parameterlist_append (self, p);
601
602 p = cpl_parameter_new_value ("gravity.preproc.remove-cosmicrays", CPL_TYPE_BOOL,
603 "Remove the cosmicrays with the statiscal method",
604 "gravity.preproc", TRUE);
605 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "remove-cosmicrays");
606 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
607 cpl_parameterlist_append (self, p);
608
609 return p;
610}
611
612cpl_parameter * gravi_parameter_add_metrology (cpl_parameterlist *self)
613{
614 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
615
616 cpl_parameter *p;
617 p = cpl_parameter_new_value ("gravity.metrology.acq-correction-delay",
618 CPL_TYPE_DOUBLE,
619 "Delay between the end of ACQ frame and correction\n "
620 "offset seen by the metrology diodes, in seconds.",
621 "gravity.metrology", 0.25);
622 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "acq-correction-delay");
623 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
624 cpl_parameterlist_append (self, p);
625
626 p = cpl_parameter_new_value ("gravity.metrology.use-fiber-dxy", CPL_TYPE_BOOL,
627 "Use the fiber position when computing OPD_TEL_CORR.",
628 "gravity.metrology", FALSE);
629 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "use-fiber-dxy");
630 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
631 cpl_parameterlist_append (self, p);
632
633 p = cpl_parameter_new_value ("gravity.metrology.use-met-rtc", CPL_TYPE_BOOL,
634 "Reduce metrology voltage with the real time algorithm\n"
635 "instead of using the pipeline’s algorithm.",
636 "gravity.metrology", FALSE);
637 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "use-met-rtc");
638 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
639 cpl_parameterlist_append (self, p);
640
641
642 p = cpl_parameter_new_value ("gravity.metrology.use-faint-met", CPL_TYPE_BOOL,
643 "In FAINT also the faint parts of the metrology \n"
644 "are used, not only the bright periods.",
645 "gravity.metrology", TRUE);
646 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "use-faint-met");
647 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
648 cpl_parameterlist_append (self, p);
649
650 p = cpl_parameter_new_value ("gravity.metrology.preswitch-delay", CPL_TYPE_INT,
651 "Delay where metrology values are ignored before\n"
652 "laser brightness is switched in faint mode, ms.",
653 "gravity.metrology", 60);
654 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "preswitch-delay");
655 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
656 cpl_parameterlist_append (self, p);
657
658 p = cpl_parameter_new_value ("gravity.metrology.postswitch-delay", CPL_TYPE_INT,
659 "Delay where metrology values are ignored after\n"
660 "laser brightness is switched in faint mode, ms.",
661 "gravity.metrology", 320);
662 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "postswitch-delay");
663 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
664 cpl_parameterlist_append (self, p);
665
666 p = cpl_parameter_new_value ("gravity.metrology.demodulate-metrology", CPL_TYPE_BOOL,
667 "Perform demodulation on the raw metrology data.",
668 "gravity.metrology", CPL_TRUE);
669 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "demodulate-metrology");
670 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
671 cpl_parameterlist_append (self, p);
672
673 p = cpl_parameter_new_value ("gravity.metrology.use-dark-offsets", CPL_TYPE_BOOL,
674 "Use diode zeros measured from the DARK when demodulating metrology.",
675 "gravity.metrology", CPL_TRUE);
676 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "use-dark-offsets");
677 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
678 cpl_parameterlist_append (self, p);
679
680 return p;
681}
682
683cpl_parameter * gravi_parameter_add_extract (cpl_parameterlist *self)
684{
685 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
686
687 cpl_parameter *p;
688 p = cpl_parameter_new_value ("gravity.preproc.ditshift-sc", CPL_TYPE_INT,
689 "Shift the time of SC DITs by an integer value to\n "
690 "account for lost frames in exposure (issue on the\n "
691 "instrument side, report to instrument team). The\n "
692 "time of all DITs in exposure are increased by\n "
693 "ditshift x PERIOD. ditshift can be 0,\n "
694 "positive (system has lost one SC DIT), or negative "
695 "(SC desynchronized).",
696 "gravity.preproc",0);
697 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "ditshift-sc");
698 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
699 cpl_parameterlist_append (self, p);
700
701 p = cpl_parameter_new_value ("gravity.preproc.extra-pixel-ft", CPL_TYPE_BOOL,
702 "Include the 6th pixels ot the FT",
703 "gravity.preproc", TRUE);
704 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "extra-pixel-ft");
705 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
706 cpl_parameterlist_append (self, p);
707
708 return p;
709}
710
711cpl_parameter * gravi_parameter_add_average_vis (cpl_parameterlist *self)
712{
713 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
714 cpl_parameter *p;
715 p = cpl_parameter_new_value ("gravity.postprocess.average-vis", CPL_TYPE_BOOL,
716 "Average the results from the different input files (if any)\n "
717 "in the output product, instead of simply appending them.",
718 "gravity.postprocess", FALSE);
719 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "average-vis");
720 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
721 cpl_parameterlist_append (self, p);
722
723 return p;
724}
725
726cpl_parameter * gravi_parameter_copy_fluxdata (cpl_parameterlist *self)
727{
728 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
729 cpl_parameter *p;
730 p = cpl_parameter_new_value ("gravity.postprocess.copy-fluxdata", CPL_TYPE_BOOL,
731 "Duplicate FLUX into FLUXDATA for OIFITS2\n "
732 "gravity.postprocess", FALSE);
733 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "copy-fluxdata");
734 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
735 cpl_parameterlist_append (self, p);
736
737 return p;
738}
739
740cpl_parameter * gravi_parameter_add_force_uncertainties (cpl_parameterlist *self)
741{
742 cpl_ensure (self, CPL_ERROR_NULL_INPUT, NULL);
743
744 cpl_parameter *p;
745 p = cpl_parameter_new_value ("gravity.postprocess.fluxerr-sc", CPL_TYPE_DOUBLE,
746 "Force the uncertainty in FLUX of SC",
747 "gravity.postprocess", 0.0);
748 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "fluxerr-sc");
749 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
750 cpl_parameterlist_append (self, p);
751
752 p = cpl_parameter_new_value ("gravity.postprocess.visamperr-sc", CPL_TYPE_DOUBLE,
753 "Force the uncertainty in VISAMP of SC",
754 "gravity.postprocess", 0.0);
755 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "visamperr-sc");
756 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
757 cpl_parameterlist_append (self, p);
758
759 p = cpl_parameter_new_value ("gravity.postprocess.visphierr-sc", CPL_TYPE_DOUBLE,
760 "Force the uncertainty in VISPHI of SC",
761 "gravity.postprocess", 0.0);
762 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "visphierr-sc");
763 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
764 cpl_parameterlist_append (self, p);
765
766 p = cpl_parameter_new_value ("gravity.postprocess.vis2err-sc", CPL_TYPE_DOUBLE,
767 "Force the uncertainty in VIS2 of SC",
768 "gravity.postprocess", 0.0);
769 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "vis2err-sc");
770 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
771 cpl_parameterlist_append (self, p);
772
773 return p;
774}
775
776cpl_error_code gravi_parameter_add_compute_snr (cpl_parameterlist *self)
777{
778 cpl_ensure_code (self, CPL_ERROR_NULL_INPUT);
779 cpl_parameter *p;
780
781 /* chi2r */
782 p = cpl_parameter_new_value ("gravity.signal.chi2r-threshold", CPL_TYPE_DOUBLE,
783 "Threshold in chi2r of the fringe-fit to declare\n "
784 "a bad value. This is usefull to detect outliers\n "
785 "or cosmic in individual frames",
786 "gravity.signal", 50.0);
787 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "chi2r-threshold");
788 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
789 cpl_parameterlist_append (self, p);
790
791 p = cpl_parameter_new_value ("gravity.signal.chi2r-sigma", CPL_TYPE_DOUBLE,
792 "Threshold in chi2r of the fringe-fit (in unit of the \n "
793 "the std of chi2r in the spectral direction) to declare\n "
794 "a bad value. This is usefull to detect outliers\n "
795 "or cosmic in individual frames",
796 "gravity.signal", 100.0);
797 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "chi2r-sigma");
798 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
799 cpl_parameterlist_append (self, p);
800
801 /* Number of FT samples to average to compute SNR and GDELAY */
802 p = cpl_parameter_new_value ("gravity.signal.nsmooth-snr-ft", CPL_TYPE_INT,
803 "Number of samples to average coherently when computing\n "
804 "the real-time SNR and GDELAY of the FT (shall correspond\n "
805 "to the atmospheric coherence time). The integration\n "
806 "window runs from -nsmooth -> +nsmooth.",
807 "gravity.signal", 5);
808 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "nsmooth-snr-ft");
809 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
810 cpl_parameterlist_append (self, p);
811
812 return CPL_ERROR_NONE;
813}
814
815cpl_error_code gravi_parameter_add_compute_signal (cpl_parameterlist *self)
816{
817 cpl_ensure_code (self, CPL_ERROR_NULL_INPUT);
818 cpl_parameter *p;
819
820 /* Maxdeg */
821 p = cpl_parameter_new_value ("gravity.signal.phase-ref-sc-maxdeg", CPL_TYPE_INT,
822 "Maximum deg for the fit of PHASE_REF",
823 "gravity.signal", 3);
824 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "phase-ref-sc-maxdeg");
825 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
826 cpl_parameterlist_append (self, p);
827
828 /* Flag to activate metrology zero calculation */
829 p = cpl_parameter_new_value ("gravity.signal.use-met-zero", CPL_TYPE_BOOL,
830 "Flag to add a constant value to OPD_DISP.\n "
831 "This constant value is taken from the header.",
832 "gravity.signal", FALSE);
833 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "use-met-zero");
834 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
835 cpl_parameterlist_append (self, p);
836
837 /* Flag to control the use of metrology in IMAGING_REF calculation */
838 p = cpl_parameter_new_enum ("gravity.signal.imaging-ref-met", CPL_TYPE_STRING,
839 "Metrology source used for IMAGING_REF calculation:\n "
840 "Use fibre coupler metrology (FC);\n "
841 "Use fibre coupler metrology corrected from pupil motion (FC_CORR);\n "
842 "Use telescope metrology (TEL).",
843 "gravity.signal", "FC", 3,
844 "FC", "FC_CORR", "TEL");
845 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "imaging-ref-met");
846 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
847 cpl_parameterlist_append (self, p);
848
849 return CPL_ERROR_NONE;
850}
851
852/*----------------------------------------------------------------------------*/
860/*----------------------------------------------------------------------------*/
861cpl_error_code gravi_parameter_add_rejection (cpl_parameterlist *self, int isCalib)
862{
863 cpl_ensure_code (self, CPL_ERROR_NULL_INPUT);
864
865 cpl_parameter *p;
866 p = cpl_parameter_new_value ("gravity.signal.snr-min-ft", CPL_TYPE_DOUBLE,
867 "SNR threshold to accept FT frames (>0). It raises the first bit (<<0)\n "
868 "of column REJECTION_FLAG of FT.",
869 "gravity.signal", isCalib ? 30.0 : 3.0);
870 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "snr-min-ft");
871 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
872 cpl_parameterlist_append (self, p);
873
874 /* OPDC_STATE threshold for fringe DET in FT */
875 p = cpl_parameter_new_value ("gravity.signal.global-state-min-ft", CPL_TYPE_DOUBLE,
876 "Minimum OPDC state to accept FT frames (>=0) It raises the second bit\n "
877 "(<<1) of column REJECTION_FLAG of FT.",
878 "gravity.signal", 2.0);
879 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "global-state-min-ft");
880 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
881 cpl_parameterlist_append (self, p);
882
883 p = cpl_parameter_new_value ("gravity.signal.global-state-max-ft", CPL_TYPE_DOUBLE,
884 "Maximum OPDC state to accept FT frames (>=0) It raises the second bit\n "
885 "(<<1) of column REJECTION_FLAG of FT.",
886 "gravity.signal", 4.0);
887 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "global-state-max-ft");
888 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
889 cpl_parameterlist_append (self, p);
890
891 /* STATE threshold for fringe DET in FT */
892 p = cpl_parameter_new_value ("gravity.signal.state-min-ft", CPL_TYPE_DOUBLE,
893 "Minimum OPDC state per baseline to accept FT frames (>=0) It raises\n "
894 "the second bit (<<1) of column REJECTION_FLAG of FT.",
895 "gravity.signal", 1.0);
896 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "state-min-ft");
897 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
898 cpl_parameterlist_append (self, p);
899
900 /* Minimum detection ratio to accept SC frame */
901 p = cpl_parameter_new_value ("gravity.signal.tracking-min-sc", CPL_TYPE_DOUBLE,
902 "Minimum ratio of accepted FT frames in order to accept a SC frames (0..1),\n "
903 "that is, for each SC DIT, the fraction of the time the\n "
904 "REJECTION_FLAG of the FT is not 0.\n "
905 "It raises the first bit (<<0) of column REJECTION_FLAG of SC",
906 "gravity.signal", 0.8);
907 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "tracking-min-sc");
908 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
909 cpl_parameterlist_append (self, p);
910
911 /* vFactor threshold to accept SC frame */
912 p = cpl_parameter_new_value ("gravity.signal.vfactor-min-sc", CPL_TYPE_DOUBLE,
913 "vFactor threshold to accept SC frame (0..1).\n ",
914 "It raises the second bit (<<1) of column REJECTION_FLAG of SC",
915 "gravity.signal", isCalib ? 0.8 : 0.1);
916 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "vfactor-min-sc");
917 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
918 cpl_parameterlist_append (self, p);
919
920 /* Threshold for OPD_PUPIL and OPD_PUPIL_STDDEV */
921 p = cpl_parameter_new_value ("gravity.signal.opd-pupil-max-sc", CPL_TYPE_DOUBLE,
922 "Maximum OPD_PUPIL (abs) to accept SC frames. It raises the third bit\n "
923 "(<<2) of column REJECTION_FLAG of SC.",
924 "gravity.signal", 9999.0);
925 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "opd-pupil-max-sc");
926 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
927 cpl_parameterlist_append (self, p);
928
929 p = cpl_parameter_new_value ("gravity.signal.opd-pupil-stddev-max-sc", CPL_TYPE_DOUBLE,
930 "Maximum OPD_PUPIL_STDDEV to accept SC frames. It\n "
931 "raises the fourth bit (<<3) of REJECTION_FLAG of SC.",
932 "gravity.signal", 2.9e-7);
933 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "opd-pupil-stddev-max-sc");
934 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
935 cpl_parameterlist_append (self, p);
936
937 return CPL_ERROR_NONE;
938}
939
940cpl_error_code gravi_parameter_add_compute_vis (cpl_parameterlist *self, int isCalib)
941{
942 cpl_ensure_code (self, CPL_ERROR_NULL_INPUT);
943 cpl_parameter *p;
944
945 /* Max-frame */
946 p = cpl_parameter_new_value ("gravity.vis.max-frame", CPL_TYPE_INT,
947 "Maximum number of frames to integrate \n "
948 "coherently into an OIFITS entry",
949 "gravity.vis", 10000);
950 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "max-frame");
951 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
952 cpl_parameterlist_append (self, p);
953
954 /* Force same time for all baselines */
955 p = cpl_parameter_new_value ("gravity.vis.force-same-time", CPL_TYPE_BOOL,
956 "Force all baseline/quantities to have\n "
957 "strictly the same TIME and MJD columns",
958 "gravity.vis", FALSE);
959 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "force-same-time");
960 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
961 cpl_parameterlist_append (self, p);
962
963 /* Debias SC VIS2 */
964 p = cpl_parameter_new_value ("gravity.vis.debias-sc", CPL_TYPE_BOOL,
965 "Subtract the V2 bias from SC",
966 "gravity.vis", TRUE);
967 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "debias-sc");
968 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
969 cpl_parameterlist_append (self, p);
970
971 /* Debias FT VIS2 */
972 p = cpl_parameter_new_value ("gravity.vis.debias-ft", CPL_TYPE_BOOL,
973 "Subtract the V2 bias from FT",
974 "gravity.vis", TRUE);
975 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "debias-ft");
976 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
977 cpl_parameterlist_append (self, p);
978
979 /* Number of bootstrap */
980 p = cpl_parameter_new_value ("gravity.vis.nboot", CPL_TYPE_INT,
981 "Number of bootstraps to compute error (1..100)",
982 "gravity.vis", isCalib ? 1 : 20);
983 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "nboot");
984 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
985 cpl_parameterlist_append (self, p);
986
987 /* Visibility correction */
988 p = cpl_parameter_new_enum ("gravity.vis.vis-correction-sc", CPL_TYPE_STRING,
989 "Correction of SC visibility from losses due to long integration,\n "
990 "using the measured visibility losses with the FT (VFACTOR\n "
991 "and/or PFACTOR) or by forcing\n "
992 "the SC visibilities to match those of the FT (FORCE). Possible\n "
993 "choices are:",
994 "gravity.vis",
995 isCalib ? "NONE" : "VFACTOR", 5, "VFACTOR", "PFACTOR",
996 "VFACTOR_PFACTOR","FORCE", "NONE");
997 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "vis-correction-sc");
998 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
999 cpl_parameterlist_append (self, p);
1000
1001 p = cpl_parameter_new_enum ("gravity.vis.vis-correction-ft", CPL_TYPE_STRING,
1002 "Correction of FT visibility from losses due to long integration,\n "
1003 "using a sliding window P_FACTOR, or its square. Choices are:",
1004 "gravity.vis", "NONE", 3, "NONE", "PFACTOR", "PFACTOR_SQUARED");
1005 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "vis-correction-ft");
1006 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1007 cpl_parameterlist_append (self, p);
1008
1009 /* Width of sliding window for FT P_FACTOR */
1010 p = cpl_parameter_new_value ("gravity.vis.pfactor-window-length", CPL_TYPE_INT,
1011 "Length of the sliding window used to calculate the FT P_FACTOR.\n "
1012 "For each FT frame, the window will run from -window-length to +window-length inclusive.",
1013 "gravity.vis", 40);
1014 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "pfactor-window-length");
1015 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1016 cpl_parameterlist_append (self, p);
1017
1018 /* Phase referencing */
1019 p = cpl_parameter_new_enum ("gravity.vis.phase-ref-sc", CPL_TYPE_STRING,
1020 "Reference phase used to integrate the SC frames.\n "
1021 "Use a self-estimate of the phase, fitted by poly. (SELF_REF)\n "
1022 "Use the FT phase only, interpolated in lbd (PHASE_REF)\n "
1023 "Use the FT+MET-SEP.UV phase (IMAGING_REF).",
1024 "gravity.vis", "AUTO", 5,
1025 "SELF_REF","PHASE_REF","IMAGING_REF","AUTO","NONE");
1026 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "phase-ref-sc");
1027 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1028 cpl_parameterlist_append (self, p);
1029
1030 /* Phase cleaning in the final OIFITS */
1031 p = cpl_parameter_new_enum ("gravity.vis.output-phase-sc", CPL_TYPE_STRING,
1032 "With DIFFERENTIAL, the mean group-delay and mean\n "
1033 "phases are removed from the output VISPHI in the\n "
1034 "final OIFITS file. With ABSOLUTE, the VISPHI is\n "
1035 "kept unmodified. With SELF_VISPHI, the internal differential\n "
1036 "phase between each spectral channel and a common \n "
1037 "reference channel is computed.\n",
1038 "gravity.vis", "AUTO", 4,
1039 "DIFFERENTIAL","ABSOLUTE","AUTO", "SELF_VISPHI");
1040 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "output-phase-sc");
1041 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1042 cpl_parameterlist_append (self, p);
1043
1044 /* string in the form [124,232] giving channel subset in case SELF_VISPHI */
1045 p = cpl_parameter_new_value ("gravity.vis.output-phase-channels", CPL_TYPE_STRING,
1046 "range (string in the form [min,max]) of channels\n "
1047 "to use a SELF_VISPHI phase reference.\n",
1048 "gravity.vis", "[0,0]");
1049 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "output-phase-channels");
1050 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1051 cpl_parameterlist_append (self, p);
1052
1053 /* Number of bootstrap */
1054 p = cpl_parameter_new_value ("gravity.vis.outlier-fraction-threshold", CPL_TYPE_DOUBLE,
1055 "Flag channels with more than this fraction of the frames\n "
1056 "affected by outliers or cosmics. These are typically detected\n "
1057 "with the thresholds options in chi2 of the fringe-fit.",
1058 "gravity.vis", 0.5);
1059 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "outlier-fraction-threshold");
1060 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1061 cpl_parameterlist_append (self, p);
1062
1063 return CPL_ERROR_NONE;
1064}
1065
1066cpl_error_code gravi_parameter_add_astrometry (cpl_parameterlist *self)
1067{
1068 cpl_ensure_code (self, CPL_ERROR_NULL_INPUT);
1069
1070 cpl_parameter *p;
1071
1072 /* just use fiber position for swaps rather than computing astrometry */
1073 p = cpl_parameter_new_value("gravity.astrometry.use-swap-fiber-pos", CPL_TYPE_BOOL,
1074 "use fiber position for swap rather than computing an astrometric solution.", "gravity.astrometry", CPL_FALSE);
1075 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "use-swap-fiber-pos");
1076 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1077 cpl_parameterlist_append (self, p);
1078
1079 /* range in RA for swap astrometry */
1080 p = cpl_parameter_new_value("gravity.astrometry.ra-lim-swap", CPL_TYPE_DOUBLE,
1081 "specify the RA range (in mas) over which to search for the astrometry of the swap. Default specifies entire field of view.", "gravity.astrometry", -1.0);
1082 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "ra-lim-swap");
1083 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1084 cpl_parameterlist_append (self, p);
1085
1086 p = cpl_parameter_new_value("gravity.astrometry.nra-swap", CPL_TYPE_INT,
1087 "number of points over the RA range for the swap.", "gravity.astrometry", 50);
1088 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "nra-swap");
1089 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1090 cpl_parameterlist_append (self, p);
1091
1092 /* range in dec for swap astrometry */
1093 p = cpl_parameter_new_value("gravity.astrometry.dec-lim-swap", CPL_TYPE_DOUBLE,
1094 "specify the dec range (in mas) over which to search for the astrometry of the swap. Default specifies entire field of view.", "gravity.astrometry", -1.0);
1095 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "dec-lim-swap");
1096 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1097 cpl_parameterlist_append (self, p);
1098
1099 p = cpl_parameter_new_value("gravity.astrometry.ndec-swap", CPL_TYPE_INT,
1100 "number of points over the dec range for the swap.", "gravity.astrometry", 50);
1101 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "ndec-swap");
1102 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1103 cpl_parameterlist_append (self, p);
1104
1105 /* average over DITs before reduction for speed */
1106 p = cpl_parameter_new_value("gravity.astrometry.average-over-dits", CPL_TYPE_BOOL,
1107 "Average over DITs before reducing astrometry for speed.", "gravity.astrometry", CPL_FALSE);
1108 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "average-over-dits");
1109 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1110 cpl_parameterlist_append (self, p);
1111
1112 p = cpl_parameter_new_value("gravity.astrometry.zoom-factor", CPL_TYPE_DOUBLE,
1113 "Factor to magnify ra/dec limits by after initial fit to find precise solution.", "gravity.astrometry", 1.0);
1114 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "zoom-factor");
1115 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1116 cpl_parameterlist_append (self, p);
1117
1118 // /* target name from FITS header */
1119 // p = cpl_parameter_new_value("gravity.astrometry.target-name", CPL_TYPE_STRING,
1120 // "specify a particular target to reduce (read from FITS header)", "gravity.astrometry", NULL);
1121
1122 // /* swap target name from FITS header */
1123 // p = cpl_parameter_new_value("gravity.astrometry.swap-target-name", CPL_TYPE_STRING,
1124 // "specify a particular target for the swap calibration in off-axis mode (read from FITS header)", "gravity.astrometry", NULL);
1125
1126 /* fiber position */
1127 // p = cpl_parameter_new_value("gravity.astrometry.fiber-pos-x", CPL_TYPE_DOUBLE,
1128 // "x position (in mas) of the fiber position to restrict the reduction to", "gravity.astrometry", 0.0);
1129 // cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "fiber-pos-x");
1130 // cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1131 // cpl_parameterlist_append (self, p);
1132
1133 // p = cpl_parameter_new_value("gravity.astrometry.fiber-pos-y", CPL_TYPE_DOUBLE,
1134 // "y position (in mas) of the fiber position to restrict the reduction to", "gravity.astrometry", 0.0);
1135 // cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "fiber-pos-y");
1136 // cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1137 // cpl_parameterlist_append (self, p);
1138
1139 /* FT flux threshold */
1140 p = cpl_parameter_new_value("gravity.astrometry.ft-mean-flux", CPL_TYPE_DOUBLE,
1141 "remove all data with FT flux below this factor times the mean", "gravity.astrometry", 0.2);
1142 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "ft-mean-flux");
1143 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1144 cpl_parameterlist_append (self, p);
1145
1146 /* Calibration strategy for computing amplitude references */
1147 const char * calib_strategies[] = {"NONE", "ALL", "SELF", "SWAP", "NEAREST"};
1148 p = cpl_parameter_new_enum_from_array("gravity.astrometry.calib-strategy", CPL_TYPE_STRING,
1149 "how to calculate the reference coherent flux\n"
1150 "\t'NONE': do not use an amplitude reference\n"
1151 "\t'ALL': use all files\n"
1152 "\t'SELF': calibrate each file individually\n"
1153 "\t'SWAP': use the swap files\n"
1154 "\t'NEAREST': use the nearest two (in time) files.",
1155 "gravity.astrometry", 0, 5, calib_strategies);
1156 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "calib-strategy");
1157 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1158 cpl_parameterlist_append (self, p);
1159
1160 /* temp for debugging */
1161 // p = cpl_parameter_new_value("gravity.astrometry.wait-for-debugger", CPL_TYPE_BOOL,
1162 // "busy wait for attaching the debugger.", "gravity.astrometry", CPL_FALSE);
1163 // cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "wait-for-debugger");
1164 // cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
1165 // cpl_parameterlist_append (self, p);
1166
1167 return CPL_ERROR_NONE;
1168}
1169
1170cpl_error_code gravi_parameter_add_image (cpl_parameterlist *self)
1171{
1172 cpl_ensure_code (self, CPL_ERROR_NULL_INPUT);
1173 cpl_parameter *p;
1174
1175 /* Fill the parameters list */
1176 /* --isotropic */
1177/* p = cpl_parameter_new_value("gravi.gravity_image.isotropic_option",
1178 CPL_TYPE_BOOL, "a flag", "gravi.gravity_image", FALSE);
1179 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "isotropic");
1180 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
1181 cpl_parameterlist_append(recipe->parameters, p);
1182*/
1183
1184 /* --pixelsize */
1185 p = cpl_parameter_new_value("gravi.gravity_image.pixelsize",
1186 CPL_TYPE_DOUBLE, "size of the pixel (milliarcseconds)",
1187 "gravi.gravity_image", 0.2);
1188 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "pixelsize");
1189 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
1190 cpl_parameterlist_append(self, p);
1191
1192 /* --dim */
1193 p = cpl_parameter_new_value("gravi.gravity_image.dim",
1194 CPL_TYPE_INT, "number of pixels per side of the image",
1195 "gravi.gravity_image", 100);
1196 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "dim");
1197 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
1198 cpl_parameterlist_append(self, p);
1199
1200 /* --regul */
1201 p = cpl_parameter_new_value("gravi.gravity_image.regul",
1202 CPL_TYPE_STRING, "name of regularization method",
1203 "gravi.gravity_image", "totvar");
1204 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "regul");
1205 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
1206 cpl_parameterlist_append(self, p);
1207
1208 /* --regul_mu */
1209 p = cpl_parameter_new_value("gravi.gravity_image.regul_mu",
1210 CPL_TYPE_DOUBLE, "global regularization weight",
1211 "gravi.gravity_image", 1E4);
1212 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "regul_mu");
1213 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
1214 cpl_parameterlist_append(self, p);
1215
1216 /* --maxeval */
1217 p = cpl_parameter_new_value("gravi.gravity_image.maxeval",
1218 CPL_TYPE_INT, "maximum number of evaluations of the objective function",
1219 "gravi.gravity_image", 2000);
1220 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "maxeval");
1221 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
1222 cpl_parameterlist_append(self, p);
1223
1224 /* --timeout */
1225 p = cpl_parameter_new_value("gravi.gravity_image.timeout",
1226 CPL_TYPE_DOUBLE, "Maximum execution time of Mira process (s)",
1227 "gravi.gravity_image", 60.);
1228 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "timeout");
1229 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
1230 cpl_parameterlist_append(self, p);
1231
1232
1233 return CPL_ERROR_NONE;
1234}
1235
1236/*----------------------------------------------------------------------------*/
1248/*----------------------------------------------------------------------------*/
1249
1250cpl_frameset * gravi_frameset_extract (cpl_frameset * frameset,
1251 const char ** frame_tags,
1252 int nb_tags)
1253{
1254 cpl_ensure (frameset, CPL_ERROR_NULL_INPUT, NULL);
1255 cpl_ensure (frame_tags, CPL_ERROR_NULL_INPUT, NULL);
1256 cpl_ensure (nb_tags>0, CPL_ERROR_ILLEGAL_INPUT, NULL);
1257
1258 int nb_frame = cpl_frameset_get_size (frameset);
1259 cpl_frameset * output_frameset = cpl_frameset_new();
1260
1261 /* Loop on frames in the frameset */
1262 for (int i = 0; i < nb_frame; i++){
1263
1264 cpl_frame * frame = cpl_frameset_get_position (frameset, i);
1265 const char * frame_tag = cpl_frame_get_tag (frame) ;
1266
1267 /* Loop on requested tags */
1268 for (int j = 0; j < nb_tags; j++) {
1269 if (strcmp(frame_tag, frame_tags[j]) == 0) {
1270 cpl_frameset_insert (output_frameset, cpl_frame_duplicate(frame));
1271 break;
1272 }
1273 }
1274
1275 } /* End loop on frames*/
1276
1277 return output_frameset;
1278}
1279
1280/*----------------------------------------------------------------------------*/
1288/*----------------------------------------------------------------------------*/
1289
1290cpl_frameset * gravi_frameset_extract_p2vm_data (cpl_frameset * frameset) {
1291 const char *tags[] = {GRAVI_P2VM_RAW};
1292 return gravi_frameset_extract (frameset, tags, 1);
1293}
1294cpl_frameset * gravi_frameset_extract_disp_data (cpl_frameset * frameset) {
1295 const char *tags[] = {GRAVI_DISP_RAW};
1296 return gravi_frameset_extract (frameset, tags, 1);
1297}
1298
1299/*----------------------------------------------------------------------------*/
1307/*----------------------------------------------------------------------------*/
1308cpl_frameset * gravi_frameset_extract_dark_data (cpl_frameset * frameset) {
1309 const char *tags[] = {GRAVI_DARK_RAW};
1310 return gravi_frameset_extract (frameset, tags, 1);
1311}
1312cpl_frameset * gravi_frameset_extract_flat_data (cpl_frameset * frameset) {
1313 const char *tags[] = {GRAVI_FLAT_RAW};
1314 return gravi_frameset_extract (frameset, tags, 1);
1315}
1316cpl_frameset * gravi_frameset_extract_diamcat_map (cpl_frameset * frameset) {
1317 const char *tags[] = {GRAVI_DIAMETER_CAT};
1318 return gravi_frameset_extract (frameset, tags, 1);
1319}
1320cpl_frameset * gravi_frameset_extract_fringe_data (cpl_frameset * frameset) {
1322 return gravi_frameset_extract (frameset, tags, 4);
1323}
1324cpl_frameset * gravi_frameset_extract_p2vmred_data (cpl_frameset * frameset) {
1327 return gravi_frameset_extract (frameset, tags, 4);
1328}
1329cpl_frameset * gravi_frameset_extract_piezotf_data (cpl_frameset * frameset) {
1330 const char *tags[] = {GRAVI_PIEZOTF_RAW};
1331 return gravi_frameset_extract (frameset, tags, 1);
1332}
1333cpl_frameset * gravi_frameset_extract_sky_data (cpl_frameset * frameset) {
1334 const char *tags[] = {GRAVI_DUAL_SKY_RAW, GRAVI_SINGLE_SKY_RAW};
1335 return gravi_frameset_extract (frameset, tags, 2);
1336}
1337cpl_frameset * gravi_frameset_extract_wave_data (cpl_frameset * frameset) {
1338 const char *tags[] = {GRAVI_WAVE_RAW};
1339 return gravi_frameset_extract (frameset, tags, 1);
1340}
1341cpl_frameset * gravi_frameset_extract_wavesc_data (cpl_frameset * frameset) {
1342 const char *tags[] = {GRAVI_WAVESC_RAW};
1343 return gravi_frameset_extract (frameset, tags, 1);
1344}
1345cpl_frameset * gravi_frameset_extract_dispvis_data (cpl_frameset * frameset) {
1346 const char *tags[] = {GRAVI_DISP_VIS};
1347 return gravi_frameset_extract (frameset, tags, 1);
1348}
1349cpl_frameset * gravi_frameset_extract_disp_map (cpl_frameset * frameset) {
1350 const char *tags[] = {GRAVI_DISP_MODEL};
1351 return gravi_frameset_extract (frameset, tags, 1);
1352}
1353cpl_frameset * gravi_frameset_extract_met_pos (cpl_frameset * frameset) {
1354 const char *tags[] = {GRAVI_DIODE_POSITION};
1355 return gravi_frameset_extract (frameset, tags, 1);
1356}
1357cpl_frameset * gravi_frameset_extract_wavelamp_map (cpl_frameset * frameset) {
1358 const char *tags[] = {GRAVI_WAVELAMP_MAP};
1359 return gravi_frameset_extract (frameset, tags, 1);
1360}
1361cpl_frameset * gravi_frameset_extract_wavelamp_data (cpl_frameset * frameset) {
1362 const char *tags[] = {GRAVI_WAVELAMP_RAW};
1363 return gravi_frameset_extract (frameset, tags, 1);
1364}
1365cpl_frameset * gravi_frameset_extract_tf_calib (cpl_frameset * frameset) {
1366 const char *tags[] = {GRAVI_TF_SINGLE_CALIB, GRAVI_TF_DUAL_CALIB,
1368 return gravi_frameset_extract (frameset, tags, 5);
1369}
1370cpl_frameset * gravi_frameset_extract_vis_calib (cpl_frameset * frameset) {
1371 const char *tags[] = {GRAVI_VIS_SINGLE_CALIB, GRAVI_VIS_DUAL_CALIB,
1373 return gravi_frameset_extract (frameset, tags, 4);
1374}
1375cpl_frameset * gravi_frameset_extract_vis_science (cpl_frameset * frameset) {
1376 const char *tags[] = {GRAVI_VIS_SINGLE_SCIENCE, GRAVI_VIS_DUAL_SCIENCE};
1377 return gravi_frameset_extract (frameset, tags, 2);
1378}
1379cpl_frameset * gravi_frameset_extract_science_data (cpl_frameset * frameset) {
1380 const char *tags[] = {GRAVI_DUAL_SCIENCE_RAW, GRAVI_SINGLE_SCIENCE_RAW};
1381 return gravi_frameset_extract (frameset, tags, 2);
1382}
1383cpl_frameset * gravi_frameset_extract_astro_target (cpl_frameset * frameset) {
1384 const char *tags[] = {GRAVI_ASTRO_TARGET};
1385 return gravi_frameset_extract (frameset, tags, 1);
1386}
1387cpl_frameset * gravi_frameset_extract_astro_swap (cpl_frameset * frameset) {
1388 const char *tags[] = {GRAVI_ASTRO_SWAP};
1389 return gravi_frameset_extract (frameset, tags, 1);
1390}
1391cpl_frameset * gravi_frameset_extract_astro_phaseref (cpl_frameset * frameset) {
1392 const char *tags[] = {GRAVI_ASTRO_CAL_PHASEREF};
1393 return gravi_frameset_extract (frameset, tags, 1);
1394}
1395cpl_frameset * gravi_frameset_extract_p2vm_map (cpl_frameset * frameset) {
1396 const char *tags[] = {GRAVI_P2VM_MAP};
1397 return gravi_frameset_extract (frameset, tags, 1);
1398}
1399cpl_frameset * gravi_frameset_extract_flat_map (cpl_frameset * frameset){
1400 const char *tags[] = {GRAVI_FLAT_MAP};
1401 return gravi_frameset_extract (frameset, tags, 1);
1402}
1403cpl_frameset * gravi_frameset_extract_dark_map (cpl_frameset * frameset) {
1404 const char *tags[] = {GRAVI_DARK_MAP};
1405 return gravi_frameset_extract (frameset, tags, 1);
1406}
1407cpl_frameset * gravi_frameset_extract_wave_map (cpl_frameset * frameset) {
1408 const char *tags[] = {GRAVI_WAVE_MAP};
1409 return gravi_frameset_extract (frameset, tags, 1);
1410}
1411cpl_frameset * gravi_frameset_extract_bad_map (cpl_frameset * frameset) {
1412 const char *tags[] = {GRAVI_BAD_MAP};
1413 return gravi_frameset_extract (frameset, tags, 1);
1414}
1415cpl_frameset * gravi_frameset_extract_biasmask_map (cpl_frameset * frameset) {
1416 const char *tags[] = {GRAVI_BIASMASK_MAP};
1417 return gravi_frameset_extract (frameset, tags, 1);
1418}
1419cpl_frameset * gravi_frameset_extract_eop_map (cpl_frameset * frameset) {
1420 const char *tags[] = {GRAVI_EOP_MAP};
1421 return gravi_frameset_extract (frameset, tags, 1);
1422}
1423cpl_frameset * gravi_frameset_extract_patch (cpl_frameset * frameset) {
1424 const char *tags[] = {GRAVI_KEY_PATCH};
1425 return gravi_frameset_extract (frameset, tags, 1);
1426}
1427cpl_frameset * gravi_frameset_extract_static_param (cpl_frameset * frameset) {
1428 const char *tags[] = {GRAVI_STATIC_PARAM};
1429 return gravi_frameset_extract (frameset, tags, 1);
1430}
1431cpl_frameset * gravi_frameset_extract_wave_param (cpl_frameset * frameset) {
1432 const char *tags[] = {GRAVI_WAVE_PARAM};
1433 return gravi_frameset_extract (frameset, tags, 1);
1434}
1435cpl_frameset * gravi_frameset_extract_pca_calib (cpl_frameset * frameset) {
1436 const char *tags[] = {GRAVI_PHASE_PCA};
1437 return gravi_frameset_extract (frameset, tags, 1);
1438}
1439
1440/*---------------------------------------------------------------------------*/
1441/*
1442 * Get the parameter from the list. Provide a default in case the parameter
1443 * is NOT in the list.
1444 */
1445/*---------------------------------------------------------------------------*/
1446
1447/*----------------------------------------------------------------------------*/
1458/*----------------------------------------------------------------------------*/
1459
1460double gravi_param_get_double_default (const cpl_parameterlist * parlist, const char * name, double def)
1461{
1462 const cpl_parameter * tmp = cpl_parameterlist_find_const(parlist, name);
1463
1464 if (tmp) {
1465 return cpl_parameter_get_double (tmp);
1466 } else {
1467 cpl_msg_info (cpl_func, "Could not find the parameter '%s':, use %f", name, def);
1468 return def;
1469 }
1470}
1471
1472int gravi_param_get_int_default (const cpl_parameterlist * parlist, const char * name, int def)
1473{
1474 const cpl_parameter * tmp = cpl_parameterlist_find_const(parlist, name);
1475
1476 if (tmp) {
1477 return cpl_parameter_get_int (tmp);
1478 } else {
1479 cpl_msg_info (cpl_func, "Could not find the parameter '%s': use %i", name, def);
1480 return def;
1481 }
1482}
1483
1484int gravi_param_get_bool_default (const cpl_parameterlist * parlist, const char * name, int def)
1485{
1486 const cpl_parameter * tmp = cpl_parameterlist_find_const(parlist, name);
1487
1488 if (tmp) {
1489 return cpl_parameter_get_bool (tmp);
1490 } else {
1491 cpl_msg_info (cpl_func, "Could not find the boolean parameter '%s': use %s", name, (def==0?"FALSE":"TRUE"));
1492 return def;
1493 }
1494}
1495
1496const char * gravi_param_get_string_default (const cpl_parameterlist * parlist, const char * name, const char * def)
1497{
1498 const cpl_parameter * tmp = cpl_parameterlist_find_const(parlist, name);
1499
1500 if (tmp) {
1501 return cpl_parameter_get_string (tmp);
1502 } else {
1503 cpl_msg_info (cpl_func, "Could not find the string parameter '%s': use %s", name, def);
1504 return def;
1505 }
1506}
1507
1508double gravi_param_get_double (const cpl_parameterlist * parlist, const char * name)
1509{
1510 const cpl_parameter * tmp = cpl_parameterlist_find_const(parlist, name);
1511 double def = 0.0;
1512
1513 if (tmp) {
1514 return cpl_parameter_get_double (tmp);
1515 } else {
1516 cpl_error_set_message (cpl_func, CPL_ERROR_ILLEGAL_INPUT, "Could not find the parameter '%s':, use %f", name, def);
1517 return def;
1518 }
1519}
1520
1521int gravi_param_get_int (const cpl_parameterlist * parlist, const char * name)
1522{
1523 const cpl_parameter * tmp = cpl_parameterlist_find_const(parlist, name);
1524 int def = 0;
1525
1526 if (tmp) {
1527 return cpl_parameter_get_int (tmp);
1528 } else {
1529 cpl_error_set_message (cpl_func, CPL_ERROR_ILLEGAL_INPUT, "Could not find the parameter '%s': use %i", name, def);
1530 return def;
1531 }
1532}
1533
1534int gravi_param_get_bool (const cpl_parameterlist * parlist, const char * name)
1535{
1536 const cpl_parameter * tmp = cpl_parameterlist_find_const(parlist, name);
1537 int def = 0;
1538
1539 if (tmp) {
1540 return cpl_parameter_get_bool (tmp);
1541 } else {
1542 cpl_error_set_message (cpl_func, CPL_ERROR_ILLEGAL_INPUT, "Could not find the boolean parameter '%s': use %s", name, (def==0?"FALSE":"TRUE"));
1543 return def;
1544 }
1545}
1546
1547const char * gravi_param_get_string (const cpl_parameterlist * parlist, const char * name)
1548{
1549 const cpl_parameter * tmp = cpl_parameterlist_find_const(parlist, name);
1550
1551 if (tmp) {
1552 return cpl_parameter_get_string (tmp);
1553 } else {
1554 cpl_error_set_message (cpl_func, CPL_ERROR_ILLEGAL_INPUT, "Could not find the string parameter '%s': use UNKNOWN", name);
1555 return "UNKNOWN";
1556 }
1557}
1558
1559cpl_error_code gravi_check_frameset (cpl_frameset *frameset, const char * tag, int min, int max)
1560{
1561 int flag = 0;
1562 int nf = cpl_frameset_count_tags (frameset, tag);
1563 char * msg = cpl_sprintf ("Need %i<#<%i '%s' in frameset (%i provided)", min, max, tag, nf);
1564
1565 if (nf < min || nf > max) {
1566 cpl_msg_error (cpl_func, "%s",msg);
1567 cpl_error_set_message (cpl_func, CPL_ERROR_ILLEGAL_INPUT, "%s", msg);
1568 flag ++;
1569 } else {
1570 cpl_msg_info (cpl_func, "%s", msg);
1571 }
1572
1573 cpl_free (msg);
1574 return (flag) ? CPL_ERROR_ILLEGAL_INPUT : CPL_ERROR_NONE;
1575}
1576
#define GRAVI_VISPHI_SINGLE_CALIB
Definition: gravi_dfs.h:115
#define GRAVI_EOP_MAP
Definition: gravi_dfs.h:78
#define GRAVI_P2VMRED_SINGLE_CALIB
Definition: gravi_dfs.h:63
#define GRAVI_TF_SINGLE_CALIB
Definition: gravi_dfs.h:106
#define GRAVI_DARK_RAW
Definition: gravi_dfs.h:46
#define GRAVI_SINGLE_SCIENCE_RAW
Definition: gravi_dfs.h:52
#define GRAVI_TF_DUAL_CALIB
Definition: gravi_dfs.h:111
#define GRAVI_PREPROC
Definition: gravi_dfs.h:60
#define GRAVI_PIEZOTF_MAP
Definition: gravi_dfs.h:86
#define GRAVI_SINGLE_SKY_RAW
Definition: gravi_dfs.h:56
#define GRAVI_ASTRO_TARGET
Definition: gravi_dfs.h:122
#define GRAVI_P2VM_MAP
Definition: gravi_dfs.h:76
#define GRAVI_NAB_CAL
Definition: gravi_dfs.h:82
#define GRAVI_P2VMRED_DUAL_CALIB
Definition: gravi_dfs.h:65
#define GRAVI_DIODE_POSITION
Definition: gravi_dfs.h:81
#define GRAVI_VIS_DUAL_SCIENCE
Definition: gravi_dfs.h:98
#define GRAVI_VIS_SINGLE_SCIENCE
Definition: gravi_dfs.h:96
#define GRAVI_KEY_PATCH
Definition: gravi_dfs.h:68
#define GRAVI_MIRA_OUTPUT_PROCATG
Definition: gravi_dfs.h:130
#define GRAVI_DISP_RAW
Definition: gravi_dfs.h:45
#define GRAVI_BAD_MAP
Definition: gravi_dfs.h:73
#define GRAVI_DUAL_CALIB_RAW
Definition: gravi_dfs.h:53
#define GRAVI_WAVE_RAW
Definition: gravi_dfs.h:47
#define GRAVI_BIASMASK_MAP
Definition: gravi_dfs.h:84
#define GRAVI_P2VMRED_SINGLE_SCIENCE
Definition: gravi_dfs.h:64
#define GRAVI_WAVE_MAP
Definition: gravi_dfs.h:75
#define GRAVI_ZP_CAL
Definition: gravi_dfs.h:128
#define GRAVI_DISP_MODEL
Definition: gravi_dfs.h:79
#define GRAVI_DUAL_SKY_RAW
Definition: gravi_dfs.h:55
#define GRAVI_FLAT_MAP
Definition: gravi_dfs.h:74
#define GRAVI_WAVE_PARAM
Definition: gravi_dfs.h:70
#define GRAVI_WAVELAMP_MAP
Definition: gravi_dfs.h:85
#define GRAVI_PIEZOTF_RAW
Definition: gravi_dfs.h:43
#define GRAVI_DIAMETER_CAT
Definition: gravi_dfs.h:80
#define GRAVI_DARK_MAP
Definition: gravi_dfs.h:77
#define GRAVI_STATIC_PARAM
Definition: gravi_dfs.h:69
#define GRAVI_P2VM_RAW
Definition: gravi_dfs.h:44
#define GRAVI_SINGLE_CALIB_RAW
Definition: gravi_dfs.h:51
#define GRAVI_VIS_DUAL_CALIB
Definition: gravi_dfs.h:99
#define GRAVI_ASTRO_SWAP
Definition: gravi_dfs.h:123
#define GRAVI_DUAL_SCIENCE_RAW
Definition: gravi_dfs.h:54
#define GRAVI_VIS_SINGLE_CALIBRATED
Definition: gravi_dfs.h:101
#define GRAVI_MIRA_INPUT_PROCATG
Definition: gravi_dfs.h:129
#define GRAVI_VIS_SINGLE_CALIB
Definition: gravi_dfs.h:97
#define GRAVI_TF_DUAL_SCIENCE
Definition: gravi_dfs.h:110
#define GRAVI_FLAT_RAW
Definition: gravi_dfs.h:50
#define GRAVI_TF_SINGLE_SCIENCE
Definition: gravi_dfs.h:107
#define GRAVI_PHASE_PCA
Definition: gravi_dfs.h:87
#define GRAVI_ASTRO_CAL_PHASEREF
Definition: gravi_dfs.h:124
#define GRAVI_WAVELAMP_RAW
Definition: gravi_dfs.h:49
#define GRAVI_P2VMRED_DUAL_SCIENCE
Definition: gravi_dfs.h:66
#define GRAVI_VISPHI_DUAL_CALIB
Definition: gravi_dfs.h:116
#define GRAVI_DISP_VIS
Definition: gravi_dfs.h:59
#define GRAVI_VISPHI_TF_DUAL_CALIB
Definition: gravi_dfs.h:118
#define GRAVI_WAVESC_RAW
Definition: gravi_dfs.h:48
#define GRAVI_VIS_DUAL_CALIBRATED
Definition: gravi_dfs.h:102
#define GRAVI_VISPHI_TF_SINGLE_CALIB
Definition: gravi_dfs.h:117
cpl_msg_info(cpl_func, "Compute WAVE_SCAN for %s", GRAVI_TYPE(type_data))
cpl_parameter * gravi_parameter_add_debug_file(cpl_parameterlist *self)
Definition: gravi_dfs.c:476
cpl_frameset * gravi_frameset_extract_wavesc_data(cpl_frameset *frameset)
Definition: gravi_dfs.c:1341
cpl_frameset * gravi_frameset_extract_met_pos(cpl_frameset *frameset)
Definition: gravi_dfs.c:1353
cpl_parameter * gravi_parameter_add_astro_file(cpl_parameterlist *self)
Definition: gravi_dfs.c:566
cpl_frameset * gravi_frameset_extract_vis_calib(cpl_frameset *frameset)
Definition: gravi_dfs.c:1370
cpl_parameter * gravi_parameter_add_profile(cpl_parameterlist *self)
Add profile parameters to the input parameter list.
Definition: gravi_dfs.c:343
int gravi_param_get_bool(const cpl_parameterlist *parlist, const char *name)
Definition: gravi_dfs.c:1534
cpl_parameter * gravi_parameter_add_static_name(cpl_parameterlist *self)
Definition: gravi_dfs.c:461
cpl_error_code gravi_parameter_add_astrometry(cpl_parameterlist *self)
Definition: gravi_dfs.c:1066
cpl_parameter * gravi_parameter_add_p2vmred_file(cpl_parameterlist *self)
Definition: gravi_dfs.c:536
cpl_parameter * gravi_parameter_add_metrology(cpl_parameterlist *self)
Definition: gravi_dfs.c:612
cpl_error_code gravi_parameter_add_rejection(cpl_parameterlist *self, int isCalib)
Add rejection parameters to the input parameter list.
Definition: gravi_dfs.c:861
cpl_frameset * gravi_frameset_extract_vis_science(cpl_frameset *frameset)
Definition: gravi_dfs.c:1375
cpl_frameset * gravi_frameset_extract_fringe_data(cpl_frameset *frameset)
Definition: gravi_dfs.c:1320
cpl_parameter * gravi_parameter_copy_fluxdata(cpl_parameterlist *self)
Definition: gravi_dfs.c:726
int gravi_param_get_bool_default(const cpl_parameterlist *parlist, const char *name, int def)
Definition: gravi_dfs.c:1484
cpl_error_code gravi_parameter_add_image(cpl_parameterlist *self)
Definition: gravi_dfs.c:1170
cpl_error_code gravi_parameter_add_compute_snr(cpl_parameterlist *self)
Definition: gravi_dfs.c:776
cpl_parameter * gravi_parameter_add_pcacalib(cpl_parameterlist *self)
Add pca calibration parameters to the input parameter list.
Definition: gravi_dfs.c:244
cpl_parameter * gravi_parameter_add_extract(cpl_parameterlist *self)
Definition: gravi_dfs.c:683
cpl_parameter * gravi_parameter_add_biasmethod(cpl_parameterlist *self)
Definition: gravi_dfs.c:581
cpl_frameset * gravi_frameset_extract_patch(cpl_frameset *frameset)
Definition: gravi_dfs.c:1423
cpl_frameset * gravi_frameset_extract_astro_target(cpl_frameset *frameset)
Definition: gravi_dfs.c:1383
cpl_parameter * gravi_parameter_add_vis_file(cpl_parameterlist *self)
Definition: gravi_dfs.c:551
cpl_parameter * gravi_parameter_add_average_vis(cpl_parameterlist *self)
Definition: gravi_dfs.c:711
cpl_frameset * gravi_frameset_extract_wave_map(cpl_frameset *frameset)
Definition: gravi_dfs.c:1407
const char * gravi_param_get_string(const cpl_parameterlist *parlist, const char *name)
Definition: gravi_dfs.c:1547
cpl_frameset * gravi_frameset_extract_astro_swap(cpl_frameset *frameset)
Definition: gravi_dfs.c:1387
cpl_frameset * gravi_frameset_extract_disp_data(cpl_frameset *frameset)
Definition: gravi_dfs.c:1294
cpl_frameset * gravi_frameset_extract_dispvis_data(cpl_frameset *frameset)
Definition: gravi_dfs.c:1345
cpl_parameter * gravi_parameter_add_wave(cpl_parameterlist *self)
Add wavelength calibration parameters to the input parameter list.
Definition: gravi_dfs.c:420
cpl_frameset * gravi_frameset_extract_bad_map(cpl_frameset *frameset)
Definition: gravi_dfs.c:1411
cpl_parameter * gravi_parameter_add_pca(cpl_parameterlist *self)
Add pca parameters to the input parameter list.
Definition: gravi_dfs.c:321
cpl_error_code gravi_parameter_disable(cpl_parameter *p)
Disable a parameter.
Definition: gravi_dfs.c:173
cpl_frameset * gravi_frameset_extract_piezotf_data(cpl_frameset *frameset)
Definition: gravi_dfs.c:1329
cpl_frameset * gravi_frameset_extract_tf_calib(cpl_frameset *frameset)
Definition: gravi_dfs.c:1365
cpl_error_code gravi_check_frameset(cpl_frameset *frameset, const char *tag, int min, int max)
Definition: gravi_dfs.c:1559
cpl_frameset * gravi_frameset_extract_science_data(cpl_frameset *frameset)
Definition: gravi_dfs.c:1379
cpl_frameset * gravi_frameset_extract_flat_data(cpl_frameset *frameset)
Definition: gravi_dfs.c:1312
cpl_parameter * gravi_parameter_add_spectrum_file(cpl_parameterlist *self)
Definition: gravi_dfs.c:506
cpl_frameset * gravi_frameset_extract(cpl_frameset *frameset, const char **frame_tags, int nb_tags)
Extract a list of tags from a frameset.
Definition: gravi_dfs.c:1250
int gravi_param_get_int(const cpl_parameterlist *parlist, const char *name)
Definition: gravi_dfs.c:1521
cpl_frameset * gravi_frameset_extract_pca_calib(cpl_frameset *frameset)
Definition: gravi_dfs.c:1435
cpl_frameset * gravi_frameset_extract_biasmask_map(cpl_frameset *frameset)
Definition: gravi_dfs.c:1415
cpl_frameset * gravi_frameset_extract_wave_data(cpl_frameset *frameset)
Definition: gravi_dfs.c:1337
cpl_parameter * gravi_parameter_add_preproc(cpl_parameterlist *self)
Add preprocessing parameters to the input parameter list.
Definition: gravi_dfs.c:394
cpl_parameter * gravi_parameter_add_force_uncertainties(cpl_parameterlist *self)
Definition: gravi_dfs.c:740
int gravi_param_get_int_default(const cpl_parameterlist *parlist, const char *name, int def)
Definition: gravi_dfs.c:1472
cpl_frameset * gravi_frameset_extract_dark_map(cpl_frameset *frameset)
Definition: gravi_dfs.c:1403
cpl_frameset * gravi_frameset_extract_p2vm_data(cpl_frameset *frameset)
Extract P2VM_RAW frame from the input frameset.
Definition: gravi_dfs.c:1290
const char * gravi_param_get_string_default(const cpl_parameterlist *parlist, const char *name, const char *def)
Definition: gravi_dfs.c:1496
cpl_frameset * gravi_frameset_extract_wavelamp_map(cpl_frameset *frameset)
Definition: gravi_dfs.c:1357
cpl_error_code gravi_dfs_set_groups(cpl_frameset *set)
Set the group as RAW or CALIB in a frameset.
Definition: gravi_dfs.c:78
cpl_parameter * gravi_parameter_add_preproc_file(cpl_parameterlist *self)
Definition: gravi_dfs.c:521
cpl_frameset * gravi_frameset_extract_p2vmred_data(cpl_frameset *frameset)
Definition: gravi_dfs.c:1324
cpl_frameset * gravi_frameset_extract_p2vm_map(cpl_frameset *frameset)
Definition: gravi_dfs.c:1395
cpl_frameset * gravi_frameset_extract_sky_data(cpl_frameset *frameset)
Definition: gravi_dfs.c:1333
cpl_parameter * gravi_parameter_add_badpix(cpl_parameterlist *self)
Add badpix parameters to the input parameter list.
Definition: gravi_dfs.c:191
cpl_frameset * gravi_frameset_extract_static_param(cpl_frameset *frameset)
Definition: gravi_dfs.c:1427
cpl_frameset * gravi_frameset_extract_wave_param(cpl_frameset *frameset)
Definition: gravi_dfs.c:1431
double gravi_param_get_double(const cpl_parameterlist *parlist, const char *name)
Definition: gravi_dfs.c:1508
cpl_error_code gravi_parameter_add_compute_vis(cpl_parameterlist *self, int isCalib)
Definition: gravi_dfs.c:940
cpl_frameset * gravi_frameset_extract_wavelamp_data(cpl_frameset *frameset)
Definition: gravi_dfs.c:1361
cpl_parameter * gravi_parameter_add_biassub_file(cpl_parameterlist *self)
Definition: gravi_dfs.c:491
cpl_frameset * gravi_frameset_extract_astro_phaseref(cpl_frameset *frameset)
Definition: gravi_dfs.c:1391
cpl_frameset * gravi_frameset_extract_flat_map(cpl_frameset *frameset)
Definition: gravi_dfs.c:1399
cpl_frameset * gravi_frameset_extract_eop_map(cpl_frameset *frameset)
Definition: gravi_dfs.c:1419
cpl_frameset * gravi_frameset_extract_diamcat_map(cpl_frameset *frameset)
Definition: gravi_dfs.c:1316
cpl_frameset * gravi_frameset_extract_disp_map(cpl_frameset *frameset)
Definition: gravi_dfs.c:1349
void gravity_print_banner(void)
Definition: gravi_dfs.c:61
cpl_frameset * gravi_frameset_extract_dark_data(cpl_frameset *frameset)
Extract DARK_RAW frame from the input frameset.
Definition: gravi_dfs.c:1308
cpl_error_code gravi_parameter_add_compute_signal(cpl_parameterlist *self)
Definition: gravi_dfs.c:815
double gravi_param_get_double_default(const cpl_parameterlist *parlist, const char *name, double def)
Get the parameter from the parameter list.
Definition: gravi_dfs.c:1460