GRAVI Pipeline Reference Manual 1.9.6
Loading...
Searching...
No Matches
gravity_vis.c
Go to the documentation of this file.
1/* $Id: gravity_vis.c,v 1.29 2011/12/3 09:16:12 nazouaoui 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
21/*
22 * $Author: nazouaoui $
23 * $Date: 2011/12/3 09:16:12 $
24 * $Revision: 1.29 $
25 * $Name: $
26 *
27 * History
28 * 12/11/2018 add static_param_frameset
29 * 26/11/2018 add static_param_frameset to call of gravi_reduce_acqcam
30 */
31
32#ifdef HAVE_CONFIG_H
33#include <config.h>
34#endif
35
36/*-----------------------------------------------------------------------------
37 Includes
38 -----------------------------------------------------------------------------*/
39
40#include <cpl.h>
41#include <stdio.h>
42#include <string.h>
43#include <time.h>
44
45#include "gravi_data.h"
46#include "gravi_pfits.h"
47#include "gravi_dfs.h"
48
49#include "gravi_utils.h"
50
51#include "gravi_calib.h"
52#include "gravi_wave.h"
53#include "gravi_p2vmred.h"
54#include "gravi_acqcam.h"
55#include "gravi_eop.h"
56#include "gravi_metrology.h"
57#include "gravi_demodulate.h"
58
59#include "gravi_signal.h"
60#include "gravi_vis.h"
61#include "gravi_tf.h"
62
63#include "gravi_preproc.h"
64
65/*-----------------------------------------------------------------------------
66 Private function prototypes
67 -----------------------------------------------------------------------------*/
68
69static int gravity_vis_create(cpl_plugin *);
70static int gravity_vis_exec(cpl_plugin *);
71static int gravity_vis_destroy(cpl_plugin *);
72static int gravity_vis(cpl_frameset *, cpl_parameterlist *);
73
74/*-----------------------------------------------------------------------------
75 Static variables
76 -----------------------------------------------------------------------------*/
77static char gravity_vis_short[] = "Compute the visibilities from raw observation of OBJECT.";
79 "This recipe is associated to the observations template. Its reduces the raw data acquired on calibrator or science targets and computes the uncalibrated visibilities, saved in an OIFITS file. If several OBJECT are provided, the recipe will reduce all of them and merge the resulting data into a single OIFITS. If several SKY_RAW are provided, the recipe reduces the first OBJECT with the first SKY file. Then each new OBJECT with the next SKY. When the number of SKYs is reached, the recipe loops back to first SKY file (so if the number of SKYs is larger than the number of OBJECTs, the last SKY won't be used). The recipe will reduce the data even if no SKY or no DARK is provided. However this will lead to wrong estimate of the visibility and squared visibility of the object. If the file DIAMETER_CAT is not provided, the recipe will use the diameter provided in the header to compute the transfer function QC parameters."
80 "\n"
81 "The tag in the DO category can be SINGLE/DUAL and CAL/SCI. They should reflect the instrument mode (SINGLE or DUAL) and the DPR.CATG of the observation (SCIENCE or CALIB). The tag in the PRO.CATG category will be SINGLE/DUAL and CAL/SCI depending on the input tag.\n"
83 "* Load the input file (loop on input OBJECT files)\n"
84 "* Extract the spectra (use BAD, DARK, SKY, FLAT files)\n"
85 "* Interpolate the spectra into a common wavelength table (use WAVE file)\n"
86 "* Compute the real-time visibilities (use P2VM file)\n"
87 "* Compute additional real-time signals (SNR, GDELAY...)\n"
88 "* Compute selection flags (= flag frames with SNR lower than threshold, vFactor lower than threshold...)\n"
89 "* Average the real-time visibilities, considering the selection flag\n"
90 "* Write the product\n"
92 GRAVI_FLAT_MAP" : flat calibration (PRO.CATG="GRAVI_FLAT_MAP")\n"
93 GRAVI_BAD_MAP" : badpixel calibration (PRO.CATG="GRAVI_BAD_MAP") \n"
94 GRAVI_WAVE_MAP" : wave calibration (PRO.CATG="GRAVI_WAVE_MAP")\n"
95 GRAVI_P2VM_MAP" : p2vm calibration (PRO.CATG="GRAVI_P2VM_MAP")\n"
96 GRAVI_DARK_MAP" : dark calibration (PRO.CATG="GRAVI_DARK_MAP")\n"
97 GRAVI_SINGLE_SCIENCE_RAW" : raw object (DPR.TYPE=OBJECT,SINGLE)\n"
98 GRAVI_SINGLE_SKY_RAW" : raw sky (DPR.TYPE=SKY,SINGLE)\n"
99 GRAVI_DISP_MODEL" (opt) : fiber dispersion model (PRO.CATG="GRAVI_DISP_MODEL")\n"
100 GRAVI_DIODE_POSITION" (opt) : met receiver position (PRO.CATG="GRAVI_DIODE_POSITION")\n"
101 GRAVI_DIAMETER_CAT" (opt) : catalog of diameter (PRO.CATG="GRAVI_DIAMETER_CAT")\n"
103 GRAVI_VIS_SINGLE_SCIENCE" : OIFITS file with uncalibrated visibilities\n"
104 GRAVI_SINGLE_SKY_MAP" (opt) : sky map\n"
105 GRAVI_P2VMRED_SINGLE_SCIENCE" (opt) : intermediate product (see detailled description of data)\n"
106 GRAVI_SPECTRUM" (opt) : intermediate product (see detailled description of data)\n"
107 GRAVI_PREPROC" (opt) : intermediate product (see detailled description of data)\n"
108 "";
109
110/*-----------------------------------------------------------------------------
111 Function code
112 -----------------------------------------------------------------------------*/
113
114/*----------------------------------------------------------------------------*/
124/*----------------------------------------------------------------------------*/
125int cpl_plugin_get_info(cpl_pluginlist * list)
126{
127 cpl_recipe * recipe = cpl_calloc(1, sizeof *recipe );
128 cpl_plugin * plugin = &recipe->interface;
129
130 if (cpl_plugin_init(plugin,
131 CPL_PLUGIN_API,
132 GRAVI_BINARY_VERSION,
133 CPL_PLUGIN_TYPE_RECIPE,
134 "gravity_vis",
137 "Nabih Azouaoui, Vincent Lapeyrere, JB. Le Bouquin",
138 PACKAGE_BUGREPORT,
143 cpl_msg_error(cpl_func, "Plugin initialization failed");
144 (void)cpl_error_set_where(cpl_func);
145 return 1;
146 }
147
148 if (cpl_pluginlist_append(list, plugin)) {
149 cpl_msg_error(cpl_func, "Error adding plugin to list");
150 (void)cpl_error_set_where(cpl_func);
151 return 1;
152 }
153
154 return 0;
155}
156
157/*----------------------------------------------------------------------------*/
165/*----------------------------------------------------------------------------*/
166static int gravity_vis_create(cpl_plugin * plugin)
167{
168 cpl_recipe * recipe;
169 cpl_parameter * p;
170
171 /* Do not create the recipe if an error code is already set */
172 if (cpl_error_get_code() != CPL_ERROR_NONE) {
173 cpl_msg_error(cpl_func, "%s():%d: An error is already set: %s",
174 cpl_func, __LINE__, cpl_error_get_where());
175 return (int)cpl_error_get_code();
176 }
177
178 if (plugin == NULL) {
179 cpl_msg_error(cpl_func, "Null plugin");
180 cpl_ensure_code(0, (int)CPL_ERROR_NULL_INPUT);
181 }
182
183 /* Verify plugin type */
184 if (cpl_plugin_get_type(plugin) != CPL_PLUGIN_TYPE_RECIPE) {
185 cpl_msg_error(cpl_func, "Plugin is not a recipe");
186 cpl_ensure_code(0, (int)CPL_ERROR_TYPE_MISMATCH);
187 }
188
189 /* Get the recipe */
190 recipe = (cpl_recipe *)plugin;
191
192 /* Create the parameters list in the cpl_recipe object */
193 recipe->parameters = cpl_parameterlist_new();
194 if (recipe->parameters == NULL) {
195 cpl_msg_error(cpl_func, "Parameter list allocation failed");
196 cpl_ensure_code(0, (int)CPL_ERROR_ILLEGAL_OUTPUT);
197 }
198
199 /* Fill the parameters list */
200 int isCalib = 0;
201
202 /* Use static names (output_procatg.fits) */
203 gravi_parameter_add_static_name (recipe->parameters);
204
205 /* Intermediate files */
206 gravi_parameter_add_biassub_file (recipe->parameters);
207 gravi_parameter_add_spectrum_file (recipe->parameters);
208 gravi_parameter_add_preproc_file (recipe->parameters);
209 gravi_parameter_add_p2vmred_file (recipe->parameters);
210 gravi_parameter_add_astro_file (recipe->parameters);
211
212 /* PCA visphi flattening */
213 gravi_parameter_add_pca (recipe->parameters);
214
215 /* Averaging */
216 gravi_parameter_add_average_vis (recipe->parameters);
217
218 /* Bias-method */
219 gravi_parameter_add_biasmethod (recipe->parameters);
220
221 /* Extraction */
222 //gravi_parameter_add_extract (recipe->parameters); // need to get the param from the p2vm
223 gravi_parameter_add_metrology (recipe->parameters);
224 // gravi_parameter_add_preproc (recipe->parameters); // now automatically handled in the preproc
225
226 /* Snr, signal, rejectio flags, vis */
227 gravi_parameter_add_compute_snr (recipe->parameters);
228 gravi_parameter_add_compute_signal (recipe->parameters);
229 gravi_parameter_add_rejection (recipe->parameters, isCalib);
230 gravi_parameter_add_compute_vis (recipe->parameters, isCalib);
231
232 /* Correct from internal transmission */
233 p = cpl_parameter_new_value ("gravity.vis.flat-flux", CPL_TYPE_BOOL,
234 "Normalize the flux (stored in OI_FLUX binary extension) with "
235 "instrument transmission recorded in the \n"
236 "input P2VM calibration map. Consequently, the flux quantity is either "
237 "the intensity level recorded \n"
238 "in the detector, thus including the instrument transmission (FALSE); "
239 "or the intensity level at the instrument entrance (TRUE).",
240 "gravity.vis", FALSE);
241 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "flat-flux");
242 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
243 cpl_parameterlist_append (recipe->parameters, p);
244
245 /* Average sky */
246 p = cpl_parameter_new_value ("gravity.preproc.average-sky", CPL_TYPE_BOOL,
247 "Average the SKYs into a master SKY. If FALSE, the recipe loops\n "
248 "over the SKY to reduce each OBJECT with a different SKY",
249 "gravity.preproc", FALSE);
250 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "average-sky");
251 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
252 cpl_parameterlist_append (recipe->parameters, p);
253
254 /* Reduce ACQ_CAM */
255 p = cpl_parameter_new_enum ("gravity.test.reduce-acq-cam", CPL_TYPE_STRING,
256 "If TRUE, reduced ACQ_CAM images. If QC, compute only\n "
257 "the QC parameters of the field part. If FALSE, ignore\n "
258 "completely the acquisition camera.",
259 "gravity.test", "FALSE",
260 3, "TRUE", "FALSE", "QC");
261 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "reduce-acq-cam");
262 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
263 cpl_parameterlist_append (recipe->parameters, p);
264
265 /* Wave color correction */
266 p = cpl_parameter_new_value ("gravity.vis.color-wave-correction", CPL_TYPE_BOOL,
267 "If TRUE, creates a new OI_WAVELENGTH_EFF with corrected wavelength",
268 "gravity.vis", FALSE);
269 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "color-wave-correction");
270 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
271 cpl_parameterlist_append (recipe->parameters, p);
272
273 /* OIFITS format */
274 p = cpl_parameter_new_value ("gravity.vis.oifits2", CPL_TYPE_BOOL,
275 "If TRUE, the output products will be fully OIFITS2 compliant. "
276 "Note that TRUE will be the default and eventually the "
277 "parameter will be removed in the future",
278 "gravity.vis", FALSE);
279 cpl_parameter_set_alias (p, CPL_PARAMETER_MODE_CLI, "oifits2");
280 cpl_parameter_disable (p, CPL_PARAMETER_MODE_ENV);
281 cpl_parameterlist_append (recipe->parameters, p);
282
283 return 0;
284}
285
286/*----------------------------------------------------------------------------*/
292/*----------------------------------------------------------------------------*/
293static int gravity_vis_exec(cpl_plugin * plugin)
294{
295
296 cpl_recipe * recipe;
297 int recipe_status;
298 cpl_errorstate initial_errorstate = cpl_errorstate_get();
299
300
301 /* Return immediately if an error code is already set */
302 if (cpl_error_get_code() != CPL_ERROR_NONE) {
303 cpl_msg_error(cpl_func, "%s():%d: An error is already set: %s",
304 cpl_func, __LINE__, cpl_error_get_where());
305 return (int)cpl_error_get_code();
306 }
307
308 if (plugin == NULL) {
309 cpl_msg_error(cpl_func, "Null plugin");
310 cpl_ensure_code(0, (int)CPL_ERROR_NULL_INPUT);
311 }
312
313 /* Verify plugin type */
314 if (cpl_plugin_get_type(plugin) != CPL_PLUGIN_TYPE_RECIPE) {
315 cpl_msg_error(cpl_func, "Plugin is not a recipe");
316 cpl_ensure_code(0, (int)CPL_ERROR_TYPE_MISMATCH);
317 }
318
319 /* Get the recipe */
320 recipe = (cpl_recipe *)plugin;
321
322 /* Verify parameter and frame lists */
323 if (recipe->parameters == NULL) {
324 cpl_msg_error(cpl_func, "Recipe invoked with NULL parameter list");
325 cpl_ensure_code(0, (int)CPL_ERROR_NULL_INPUT);
326 }
327 if (recipe->frames == NULL) {
328 cpl_msg_error(cpl_func, "Recipe invoked with NULL frame set");
329 cpl_ensure_code(0, (int)CPL_ERROR_NULL_INPUT);
330 }
331
332 /* Invoke the recipe */
333 recipe_status = gravity_vis(recipe->frames, recipe->parameters);
334
335 /* Ensure DFS-compliance of the products */
336 if (cpl_dfs_update_product_header(recipe->frames)) {
337 if (!recipe_status){
338 recipe_status = (int)cpl_error_get_code();
339 }
340 }
341
342 if (!cpl_errorstate_is_equal(initial_errorstate)) {
343 /* Dump the error history since recipe execution start.
344 At this point the recipe cannot recover from the error */
345 cpl_errorstate_dump(initial_errorstate, CPL_FALSE, NULL);
346 }
347
348 return recipe_status;
349}
350
351/*----------------------------------------------------------------------------*/
357/*----------------------------------------------------------------------------*/
358static int gravity_vis_destroy(cpl_plugin * plugin)
359{
360 cpl_recipe * recipe;
361
362 if (plugin == NULL) {
363 cpl_msg_error(cpl_func, "Null plugin");
364 cpl_ensure_code(0, (int)CPL_ERROR_NULL_INPUT);
365 }
366
367 /* Verify plugin type */
368 if (cpl_plugin_get_type(plugin) != CPL_PLUGIN_TYPE_RECIPE) {
369 cpl_msg_error(cpl_func, "Plugin is not a recipe");
370 cpl_ensure_code(0, (int)CPL_ERROR_TYPE_MISMATCH);
371 }
372
373 /* Get the recipe */
374 recipe = (cpl_recipe *)plugin;
375
376 cpl_parameterlist_delete(recipe->parameters);
377
378 return 0;
379}
380
381
382/*----------------------------------------------------------------------------*/
390/*----------------------------------------------------------------------------*/
391static int gravity_vis(cpl_frameset * frameset,
392 cpl_parameterlist * parlist)
393{
394 cpl_frameset * recipe_frameset=NULL, * wavecalib_frameset=NULL, * dark_frameset=NULL,
395 * darkcalib_frameset=NULL, * sky_frameset=NULL, * flatcalib_frameset=NULL, * p2vmcalib_frameset=NULL,
396 * badcalib_frameset=NULL, *used_frameset=NULL, * current_frameset=NULL, * dispcalib_frameset=NULL,
397 * metpos_frameset=NULL, * diamcat_frameset = NULL, *eop_frameset = NULL, *patch_frameset = NULL,
398 * static_param_frameset=NULL, * pcacalib_frameset = NULL;
399
400 cpl_frame * frame=NULL;
401
402 const char * frame_tag=NULL;
403 char * proCatg = NULL, * mode=NULL, * redCatg = NULL, * skyCatg = NULL;
404
405 gravi_data * p2vm_map=NULL, * data=NULL, * wave_map=NULL, * dark_map=NULL,
406 * profile_map=NULL, * badpix_map=NULL, * preproc_data=NULL, * p2vmred_data=NULL, * tmpvis_data=NULL,
407 * vis_data=NULL, * disp_map=NULL, * diodepos_data=NULL, * diamcat_data=NULL, *eop_map=NULL,
408 * static_param_data=NULL, * pca_calib_data=NULL;
409 gravi_data ** sky_maps = NULL;
410 cpl_propertylist ** p2vm_qcs = NULL;
411
412 int nb_frame, nb_sky;
413 char * input_data_type;
414
415 /* Message */
418
419 cpl_ensure_code(gravi_dfs_set_groups(frameset) == CPL_ERROR_NONE,
420 cpl_error_get_code()) ;
421
422 /* Dispatch the frameset */
423 p2vmcalib_frameset = gravi_frameset_extract_p2vm_map (frameset);
424 darkcalib_frameset = gravi_frameset_extract_dark_map (frameset);
425 wavecalib_frameset = gravi_frameset_extract_wave_map (frameset);
426 dark_frameset = gravi_frameset_extract_dark_data (frameset);
427 flatcalib_frameset = gravi_frameset_extract_flat_map (frameset);
428 badcalib_frameset = gravi_frameset_extract_bad_map (frameset);
429 dispcalib_frameset = gravi_frameset_extract_disp_map (frameset);
430 pcacalib_frameset = gravi_frameset_extract_pca_calib (frameset);
431 metpos_frameset = gravi_frameset_extract_met_pos (frameset);
432 diamcat_frameset = gravi_frameset_extract_diamcat_map (frameset);
433 eop_frameset = gravi_frameset_extract_eop_map (frameset);
434 patch_frameset = gravi_frameset_extract_patch (frameset);
435 static_param_frameset = gravi_frameset_extract_static_param (frameset);
436
437 recipe_frameset = gravi_frameset_extract_fringe_data (frameset);
438 sky_frameset = gravi_frameset_extract_sky_data (frameset);
439
440 /* To use this recipe the frameset must contain the p2vm, wave and
441 * gain calibration file. */
442 if ( cpl_frameset_get_size (p2vmcalib_frameset) !=1 ||
443 cpl_frameset_get_size (wavecalib_frameset) !=1 ||
444 cpl_frameset_get_size (flatcalib_frameset) !=1 ||
445 cpl_frameset_get_size (badcalib_frameset) != 1 ||
446 cpl_frameset_get_size (recipe_frameset) < 1 ||
447 (cpl_frameset_is_empty (dark_frameset) &&
448 cpl_frameset_is_empty (darkcalib_frameset) &&
449 cpl_frameset_is_empty (sky_frameset)) ) {
450 cpl_error_set_message (cpl_func, CPL_ERROR_ILLEGAL_INPUT,
451 "Illegal number of P2VM, FLAT, WAVE, BAD, DARK or SKY, OBJECT file on the frameset");
452 cpl_error_set_message (cpl_func, CPL_ERROR_ILLEGAL_INPUT,
453 "See online help: esorex --man gravity_vis");
454 goto cleanup;
455 }
456
457 /* Force some options if phase flattening is to be performed */
458 if (gravi_param_get_bool (parlist, "gravity.vis.flatten-visphi")) {
459 cpl_parameter *phase_ref = cpl_parameterlist_find (parlist, "gravity.vis.phase-ref-sc");
460 cpl_parameter *output_phase = cpl_parameterlist_find (parlist, "gravity.vis.output-phase-sc");
461
462 if (strcmp (cpl_parameter_get_string(phase_ref), "SELF_REF") != 0) {
463 cpl_msg_warning (cpl_func, "VISPHI flattening requires phase-ref-sc=SELF_REF, forcing");
464 cpl_parameter_set_string (phase_ref, "SELF_REF");
465 }
466
467 if (strcmp (cpl_parameter_get_string(output_phase), "SELF_VISPHI") != 0) {
468 cpl_msg_warning (cpl_func, "VISPHI flattening requires output-phase-sc=SELF_VISPHI, forcing");
469 cpl_parameter_set_string (output_phase, "SELF_VISPHI");
470 }
471 }
472
473 /* Insert calibration frame into the used frameset */
474 used_frameset = cpl_frameset_new();
475
476 /*
477 * Identify the DARK in the input frameset
478 */
479
480 if (!cpl_frameset_is_empty (dark_frameset)) {
481
482 frame = cpl_frameset_get_position (dark_frameset, 0);
483 data = gravi_data_load_rawframe (frame, used_frameset);
484 gravi_data_patch (data, patch_frameset);
485 gravi_data_detector_cleanup (data, parlist);
486
487 /* Compute dark */
488 dark_map = gravi_compute_dark (data);
489 FREE (gravi_data_delete, data);
490
491 CPLCHECK_CLEAN ("Could not compute the DARK map");
492
493 /* Save the dark map */
494 gravi_data_save_new (dark_map, frameset, NULL, NULL, parlist,
495 NULL, frame, "gravity_vis",
496 NULL, GRAVI_DARK_MAP);
497
498 CPLCHECK_CLEAN ("Could not save the DARK map");
499 }
500 else if (!cpl_frameset_is_empty (darkcalib_frameset)) {
501
502 frame = cpl_frameset_get_position (darkcalib_frameset, 0);
503 dark_map = gravi_data_load_frame (frame, used_frameset);
504
505 CPLCHECK_CLEAN ("Could not load the DARK map");
506 }
507 else
508 cpl_msg_info (cpl_func, "There is no DARK in the frame set");
509
510 /* Identify the BAD in the input frameset */
511 frame = cpl_frameset_get_position (badcalib_frameset, 0);
512 badpix_map = gravi_data_load_frame (frame, used_frameset);
513
514 /* Identify the FLAT in the input frameset */
515 frame = cpl_frameset_get_position (flatcalib_frameset, 0);
516 profile_map = gravi_data_load_frame (frame, used_frameset);
517
518 /* Identify the WAVE in the input frameset */
519 frame = cpl_frameset_get_position (wavecalib_frameset, 0);
520 wave_map = gravi_data_load_frame (frame, used_frameset);
521
522 /* Identify the P2VM in the input frameset */
523 frame = cpl_frameset_get_position (p2vmcalib_frameset, 0);
524 p2vm_map = gravi_data_load_frame (frame, used_frameset);
525 cpl_parameter * param_extrapixel = gravi_pfits_get_extrapixel_param(gravi_data_get_header(p2vm_map));
526 cpl_parameterlist_append(parlist, param_extrapixel);
527
528 /* Load the DISP_MODEL in the input frameset */
529 if (!cpl_frameset_is_empty (dispcalib_frameset)) {
530 frame = cpl_frameset_get_position (dispcalib_frameset, 0);
531 disp_map = gravi_data_load_frame (frame, used_frameset);
532 }
533 else
534 cpl_msg_info (cpl_func, "There is no DISP_MODEL in the frameset");
535
536 /* Load the DIODE_POSITION in the input frameset */
537 if (!cpl_frameset_is_empty (metpos_frameset)) {
538 frame = cpl_frameset_get_position (metpos_frameset, 0);
539 diodepos_data = gravi_data_load_frame (frame, used_frameset);
540 }
541 else
542 cpl_msg_info (cpl_func, "There is no DIODE_POSITION in the frameset");
543
544 /* Load the EOP_PARAM */
545 if ( !cpl_frameset_is_empty (eop_frameset) ) {
546 frame = cpl_frameset_get_position (eop_frameset, 0);
547 eop_map = gravi_data_load_frame (frame, used_frameset);
548 }
549 else
550 cpl_msg_info (cpl_func, "There is no EOP_PARAM in the frameset");
551
552 /* START EKW 12/11/2018 read constant parameter from calibration file */
553 /* Load the STATIC_PARAM Parameter */
554 if (!cpl_frameset_is_empty (static_param_frameset)) {
555 frame = cpl_frameset_get_position (static_param_frameset, 0);
556 static_param_data = gravi_data_load_frame (frame, used_frameset);
557 }
558 else
559 cpl_msg_info (cpl_func, "There is no STATIC_PARAM in the frameset");
560
561 /* Load the DIAMETER_CAT */
562 if ( !cpl_frameset_is_empty (diamcat_frameset) ) {
563 frame = cpl_frameset_get_position (diamcat_frameset, 0);
564 diamcat_data = gravi_data_load_frame (frame, used_frameset);
565 }
566 else
567 cpl_msg_info (cpl_func, "There is no DIAMETER_CAT in the frameset");
568
569 if ( !cpl_frameset_is_empty (pcacalib_frameset)) {
570 frame = cpl_frameset_get_position (pcacalib_frameset, 0);
571 pca_calib_data = gravi_data_load_frame (frame, used_frameset);
572 } else
573 cpl_msg_info (cpl_func, "There is no PHASE_PCA in the frameset");
574
575 CPLCHECK_CLEAN ("Error while loading the calibration maps");
576
577 /*
578 * Select the PRO CATG (based on first frame)
579 */
580
581 frame_tag = cpl_frame_get_tag (cpl_frameset_get_position (recipe_frameset, 0));
582
583 if ((strcmp(frame_tag, GRAVI_DUAL_CALIB_RAW) == 0)) {
584 redCatg = cpl_sprintf (GRAVI_P2VMRED_DUAL_CALIB);
585 proCatg = cpl_sprintf (GRAVI_VIS_DUAL_CALIB);
586 skyCatg = cpl_sprintf (GRAVI_DUAL_SKY_MAP);
587 mode = cpl_sprintf ("gravi_dual");
588 input_data_type = cpl_sprintf ("raw_calibrator");
589 }
590 else if ((strcmp(frame_tag, GRAVI_DUAL_SCIENCE_RAW) == 0)) {
591 redCatg = cpl_sprintf (GRAVI_P2VMRED_DUAL_SCIENCE);
592 proCatg = cpl_sprintf (GRAVI_VIS_DUAL_SCIENCE);
593 skyCatg = cpl_sprintf (GRAVI_DUAL_SKY_MAP);
594 mode = cpl_sprintf ("gravi_dual");
595 input_data_type = cpl_sprintf ("raw_science");
596 }
597 else if ((strcmp(frame_tag, GRAVI_SINGLE_CALIB_RAW) == 0)) {
598 redCatg = cpl_sprintf (GRAVI_P2VMRED_SINGLE_CALIB);
599 proCatg = cpl_sprintf (GRAVI_VIS_SINGLE_CALIB);
600 skyCatg = cpl_sprintf (GRAVI_SINGLE_SKY_MAP);
601 mode = cpl_sprintf ("gravi_single");
602 input_data_type = cpl_sprintf ("raw_calibrator");
603 }
604 else if ((strcmp(frame_tag, GRAVI_SINGLE_SCIENCE_RAW) == 0)) {
605 redCatg = cpl_sprintf (GRAVI_P2VMRED_SINGLE_SCIENCE);
606 proCatg = cpl_sprintf (GRAVI_VIS_SINGLE_SCIENCE);
607 skyCatg = cpl_sprintf (GRAVI_SINGLE_SKY_MAP);
608 mode = cpl_sprintf ("gravi_single");
609 input_data_type = cpl_sprintf ("raw_science");
610 }
611 else {
612 cpl_error_set_message (cpl_func, CPL_ERROR_ILLEGAL_INPUT,
613 "Cannot recognize the input DO.CATG");
614 goto cleanup;
615 }
616
617 cpl_msg_info (cpl_func,"Mode of the first frame is: %s (will be used for all frames)", mode);
618
619 /*
620 * Mode for the SKY
621 */
622 int averageSky = gravi_param_get_bool (parlist,"gravity.preproc.average-sky");
623
624 /*
625 * Loop on input SKY frames to be reduced
626 */
627 nb_sky = cpl_frameset_get_size (sky_frameset);
628 sky_maps = cpl_calloc (CPL_MAX(nb_sky,1), sizeof(gravi_data*));
629
630 for (int isky = 0; isky < nb_sky; isky++){
631
632 /* Load the raw SKY */
633 cpl_msg_info (cpl_func, " ***** SKY %d over %d ***** ", isky+1, nb_sky);
634 frame = cpl_frameset_get_position (sky_frameset, isky);
635 data = gravi_data_load_rawframe (frame, used_frameset);
636 gravi_data_patch (data, patch_frameset);
637 gravi_data_detector_cleanup (data, parlist);
638
639 /* Compute the SKY map */
640 sky_maps[isky] = gravi_compute_dark (data);
641 FREE (gravi_data_delete, data);
642
643 CPLCHECK_CLEAN ("Error while computing the sky_map");
644
645 /* Save the SKY map */
646 if (averageSky == 0) {
647 char filename_suffix[20];
648 snprintf(filename_suffix, 16, "%d", isky);
649 gravi_data_save_new (sky_maps[isky], frameset, NULL, filename_suffix,
650 parlist, NULL, frame, "gravity_vis",
651 NULL, skyCatg);
652 CPLCHECK_CLEAN ("Could not save the sky");
653 }
654 }
655
656 /*
657 * Average the sky if requested
658 */
659
660 if (averageSky == 1) {
661 cpl_msg_info (cpl_func, "Do a MASTER SKY from the %d skys", nb_sky);
662
663 gravi_data * msky_map;
664 msky_map = gravi_average_dark (sky_maps, nb_sky);
665 CPLCHECK_CLEAN ("Cannot do master sky");
666
667 gravi_data_save_new (msky_map, frameset, NULL, NULL,
668 parlist, sky_frameset,
669 cpl_frameset_get_position (sky_frameset, 0),
670 "gravity_vis", NULL, skyCatg);
671 CPLCHECK_CLEAN ("Cannot save master sky");
672
673 /* Add all sky to used_frameset, and move pointers */
674 cpl_frameset_join (used_frameset, sky_frameset);
675 for (int isky = 0; isky < nb_sky; isky++)
676 FREE (gravi_data_delete, sky_maps[isky]);
677 sky_maps[0] = msky_map;
678 nb_sky = 1;
679 }
680
681 /*
682 * Loop on input RAW frames to be reduced
683 */
684
685 nb_frame = cpl_frameset_get_size (recipe_frameset);
686 p2vm_qcs = cpl_malloc(sizeof(cpl_propertylist*) * nb_frame);
687
688 for (int ivis = 0; ivis < nb_frame; ivis++){
689 p2vm_qcs[ivis] = cpl_propertylist_new();
690 int isky;
691 char filename_suffix[20];
692 snprintf(filename_suffix, 16, "%d", ivis);
693 current_frameset = cpl_frameset_duplicate (used_frameset);
694
695 cpl_msg_info (cpl_func, " ***** OBJECT %d over %d ***** ", ivis+1, nb_frame);
696
697 /*
698 * Identify the SKY for this OBJECT
699 */
700 isky = nb_sky>0 ? ivis % nb_sky : 0;
701
702 if (nb_sky == 0) {
703 /* No SKY */
704 cpl_msg_info (cpl_func, "There is no SKY in the frameset");
705 }
706 else if (averageSky) {
707 /* Use master SKY already computed, already in frameset */
708 cpl_msg_info (cpl_func, "Use MASTER SKY (already reduced)");
709 }
710 else {
711 /* SKY already computed, add in the used_frameset */
712 cpl_msg_info (cpl_func, "Use SKY %i over %i (already reduced)", isky+1, nb_sky);
713 frame = cpl_frameset_get_position (sky_frameset, isky);
714
715 /* Add this frame to the current_frameset as well */
716 cpl_frameset_insert (current_frameset, cpl_frame_duplicate (frame));
717 }
718
719 /*
720 * Reduce the OBJECT
721 */
722
723 frame = cpl_frameset_get_position (recipe_frameset, ivis);
724 /* Add this frame to the used frameset */
725 cpl_frameset_insert(used_frameset, cpl_frame_duplicate (frame));
726
727 /* Start processing */
728 data = gravi_data_load_rawframe (frame, current_frameset);
729 gravi_data_patch (data, patch_frameset);
730 gravi_data_detector_cleanup (data, parlist);
731
732 /* Option save the bias-subtracted file */
733 if (gravi_param_get_bool (parlist,"gravity.dfs.bias-subtracted-file")) {
734
735 gravi_data_save_new (data, frameset, NULL, filename_suffix, parlist,
736 current_frameset, frame, "gravity_vis",
737 NULL, "BIAS_SUBTRACTED");
738
739 CPLCHECK_CLEAN ("Cannot save the BIAS_SUBTRACTED product");
740 }
741
742 /* Check the shutters */
743 if ( !gravi_data_check_shutter_open (data) ) {
744 cpl_msg_warning (cpl_func, "Shutter problem in the OBJECT");
745 }
746
747 /* Extract spectrum */
748 preproc_data = gravi_extract_spectrum (data, profile_map, dark_map,
749 badpix_map, sky_maps[isky],
750 parlist, GRAVI_DET_ALL);
751 CPLCHECK_CLEAN ("Cannot extract spectrum");
752
753 /* Option save the spectrum file */
754 if (gravi_param_get_bool (parlist,"gravity.dfs.spectrum-file")) {
755 gravi_data_save_new (preproc_data, frameset, NULL, filename_suffix,
756 parlist, current_frameset, frame,
757 "gravity_vis", NULL, GRAVI_SPECTRUM);
758 CPLCHECK_CLEAN ("Cannot save the SPECTRUM product");
759 }
760
761 /* Rescale to common wavelength */
762 gravi_align_spectrum (preproc_data, wave_map, p2vm_map, GRAVI_DET_ALL);
763 CPLCHECK_CLEAN ("Cannot re-interpolate spectrum");
764
765 /* Option save the spectrum-aligned file */
766 if (gravi_param_get_bool (parlist,"gravity.dfs.spectrum-file")) {
767 gravi_data_save_new (preproc_data, frameset, NULL, filename_suffix,
768 parlist, current_frameset, frame,
769 "gravity_vis", NULL, GRAVI_SPECTRUM_ALIGNED);
770 CPLCHECK_CLEAN ("Cannot save the SPECTRUM_ALIGNED product");
771 }
772
773 /* Preproc the Acquisition Camera */
774 if (!strcmp (gravi_param_get_string (parlist, "gravity.test.reduce-acq-cam"), "TRUE") ||
775 !strcmp (gravi_param_get_string (parlist, "gravity.test.reduce-acq-cam"), "QC")) {
776 gravi_preproc_acqcam (preproc_data, data, badpix_map);
777 CPLCHECK_CLEAN ("Cannot preproc ACQ");
778 }
779
780 /* Option save the preproc file */
781 if (gravi_param_get_bool (parlist,"gravity.dfs.preproc-file")) {
782 gravi_data_save_new (preproc_data, frameset, NULL, filename_suffix,
783 parlist, current_frameset, frame,
784 "gravity_vis", NULL, GRAVI_PREPROC);
785 CPLCHECK_CLEAN ("Cannot save the PREPROC product");
786 }
787
788 /* Copy metrology and subtract background to preproc */
789 gravi_data_move_ext (preproc_data, data, GRAVI_METROLOGY_EXT);
790 cpl_boolean subtract_met_dark = dark_map != NULL && cpl_parameter_get_bool(
791 cpl_parameterlist_find_const(parlist, "gravity.metrology.use-dark-offsets"));
792
793 if (subtract_met_dark)
794 gravi_subtract_met_dark (preproc_data, dark_map);
795
796 /* Demodulate the metrology if requested */
797 if (gravi_param_get_bool (parlist, "gravity.metrology.demodulate-metrology")) {
798 gravi_metrology_demodulate(preproc_data, subtract_met_dark);
799 CPLCHECK_CLEAN ("Cannot demodulate metrology");
800 }
801
802 /* Move extensions from raw_data and delete it */
803 gravi_data_move_ext (preproc_data, data, GRAVI_ARRAY_GEOMETRY_EXT);
804 gravi_data_move_ext (preproc_data, data, GRAVI_OPTICAL_TRAIN_EXT);
805 gravi_data_move_ext (preproc_data, data, GRAVI_OPDC_EXT);
806 gravi_data_move_ext (preproc_data, data, GRAVI_FDDL_EXT);
807 FREE (gravi_data_delete, data);
808 CPLCHECK_CLEAN ("Cannot move ext");
809
810 /* Compute the flux and visibilities for each telescope and
811 * per acquisition with the P2VM applied to preproc_data */
812 p2vmred_data = gravi_compute_p2vmred(preproc_data, p2vm_map, mode,
813 parlist);
814 CPLCHECK_CLEAN ("Cannot apply p2vm to the preproc data");
815
816 /* Reduce the Acquisition Camera */
817 if (!strcmp (gravi_param_get_string (parlist, "gravity.test.reduce-acq-cam"), "TRUE")) {
818 int saveAcqTable = 1;
819 gravi_reduce_acqcam (p2vmred_data, preproc_data, sky_maps[isky], dark_map, static_param_data, saveAcqTable);
820 }
821 else if (!strcmp (gravi_param_get_string (parlist, "gravity.test.reduce-acq-cam"), "QC")) {
822 int saveAcqTable = 0;
823 gravi_reduce_acqcam (p2vmred_data, preproc_data, sky_maps[isky], dark_map, static_param_data, saveAcqTable);
824 }
825
826 /* Move extensions and delete preproc */
827 gravi_data_move_ext (p2vmred_data, preproc_data, GRAVI_IMAGING_DATA_ACQ_EXT);
828 gravi_data_move_ext (p2vmred_data, preproc_data, GRAVI_METROLOGY_EXT);
829 gravi_data_move_ext (p2vmred_data, preproc_data, GRAVI_FDDL_EXT);
830 gravi_data_move_ext (p2vmred_data, preproc_data, GRAVI_OPDC_EXT);
831 FREE (gravi_data_delete, preproc_data);
832 CPLCHECK_CLEAN ("Cannot delete preproc");
833
834 /* Reduce the OPDC table */
835 gravi_compute_opdc_state (p2vmred_data);
836 CPLCHECK_CLEAN ("Cannot reduce OPDC");
837
838 /* Reduce the metrology into OI_VIS_MET */
839 gravi_metrology_reduce (p2vmred_data, eop_map, static_param_data, diodepos_data, parlist);
840 CPLCHECK_CLEAN ("Cannot reduce metrology");
841
842 /* Compute the uv and pointing directions with ERFA */
843 gravi_compute_pointing_uv (p2vmred_data, eop_map);
844 CPLCHECK_CLEAN ("Cannot compute pointing");
845
846 /* Compute the QC0 about tau0 from piezo signals */
847 gravi_compute_tau0 (p2vmred_data);
848 CPLCHECK_CLEAN ("Cannot compute QC for tau0");
849
850 /* Compute QC for the Fringe Tracker injection */
851 gravi_compute_qc_injection (p2vmred_data);
852 CPLCHECK_CLEAN ("Cannot compute QC for FT injection");
853
854 /* Compute QC for the Fringe Tracker OPD calculation */
856 CPLCHECK_CLEAN ("Cannot compute QC for FT OPD estimator");
857
858 /* Find outliers */
859 gravi_compute_outliers (p2vmred_data, parlist);
860 CPLCHECK_MSG ("Cannot compute outliers");
861
862 /* Compute the SNR_BOOT and GDELAY_BOOT */
863 gravi_compute_snr (p2vmred_data, parlist);
864 CPLCHECK_MSG ("Cannot compute SNR");
865
866 /* Compute the signals for averaging */
867 gravi_compute_signals (p2vmred_data, disp_map, parlist);
868 CPLCHECK_MSG ("Cannot compute signals");
869
870 /* Temporary copy to allow averaging over all frames */
871 gravi_copy_p2vm_qcs(p2vmred_data, p2vm_qcs[ivis]);
872 CPLCHECK_MSG ("Cannot copy QC for averaging");
873
874 /* Compute rejection flags for averaging */
875 gravi_compute_rejection (p2vmred_data, parlist);
876 CPLCHECK_MSG ("Cannot compute rejection flags signals");
877
878 /* Save the p2vmreduced file */
879 if (gravi_param_get_bool (parlist,"gravity.dfs.p2vmred-file")) {
880
881 gravi_data_save_new (p2vmred_data, frameset, NULL, filename_suffix,
882 parlist, current_frameset, frame,
883 "gravity_vis", NULL, redCatg);
884
885 CPLCHECK_CLEAN ("Cannot save the P2VMREDUCED product");
886 }
887
888 /* Loop on the wanted sub-integration */
889 cpl_size current_frame = 0;
890 while (current_frame >= 0)
891 {
892 /* Visibility and flux are averaged and the followings
893 * are saved in tables VIS, VIS2 and T3 */
894 tmpvis_data = gravi_compute_vis (p2vmred_data, parlist, &current_frame);
895 CPLCHECK_CLEAN ("Cannot average the P2VMRED frames into VIS");
896
897 /* Set the mean TIME and mean MJD if required */
898 if (gravi_param_get_bool (parlist, "gravity.vis.force-same-time") ) {
899 cpl_msg_info (cpl_func,"Force same time for all quantities/baselines");
900 gravi_vis_force_time (tmpvis_data);
901 CPLCHECK_CLEAN ("Cannot average the TIME in OI_VIS");
902 }
903
904 /* Copy the acquisition camera if requested. */
905 if (current_frame < 0 && !strcmp (gravi_param_get_string (parlist, "gravity.test.reduce-acq-cam"), "TRUE"))
906 {
907 cpl_msg_info (cpl_func, "Copy ACQ into the VIS file");
909 }
910
911 /* Merge with already existing */
912 if (vis_data == NULL) {
913 vis_data = tmpvis_data; tmpvis_data = NULL;
914 }
915 else {
916 cpl_msg_info (cpl_func,"Merge with previous OI_VIS");
917 gravi_data_append (vis_data, tmpvis_data, 1);
918 FREE (gravi_data_delete, tmpvis_data);
919 }
920 }
921
922 /* Save the astro file, which is a lighter version of the p2vmreduced */
923 if (gravi_param_get_bool (parlist,"gravity.dfs.astro-file")) {
924
925 gravi_data_clean_for_astro (p2vmred_data);
926 gravi_data_save_new (p2vmred_data, frameset, NULL, filename_suffix,
927 parlist, current_frameset, frame,
928 "gravity_vis", NULL, GRAVI_ASTROREDUCED);
929
930 CPLCHECK_CLEAN ("Cannot save the ASTROREDUCED product");
931 }
932
933 cpl_msg_info (cpl_func,"Free the p2vmreduced");
934 FREE (cpl_frameset_delete, current_frameset);
935 FREE (gravi_data_delete, p2vmred_data);
936
937 }
938 /* End loop on the input files to reduce */
939
940 /* Use the PCA calibration to flatten the VISPHI */
941 if (gravi_param_get_bool (parlist, "gravity.vis.flatten-visphi")) {
942 cpl_msg_info (cpl_func, "Flatten VISPHI using PCA");
943 gravi_flatten_vis(vis_data, pca_calib_data);
944 CPLCHECK_CLEAN ("Cannot apply the VISPHI flattening");
945 }
946
947 /* Compute QC parameters */
948 cpl_msg_info (cpl_func, "Computing QC parameters for visibilities");
949 gravi_compute_vis_qc (vis_data, frameset, p2vm_qcs, nb_frame, input_data_type);
950 CPLCHECK_CLEAN ("Cannot compute VIS QCs");
951
952 /* Compute the QC parameters of the TF
953 * FIXME: compute QC TF only for CALIB star */
954 gravi_compute_tf_qc (vis_data, diamcat_data);
955
956 /* Eventually flatten the OI_FLUX */
957 if (gravi_param_get_bool (parlist, "gravity.vis.flat-flux")) {
958
959 cpl_msg_info (cpl_func, "Flatten the FLUX with the internal P2VM spectrum");
960 gravi_flat_flux (vis_data, p2vm_map);
961 CPLCHECK_CLEAN ("Cannot flat the OI_FLUX");
962
963 } else {
964 cpl_msg_info (cpl_func, "Don't flatten the FLUX with the internal P2VM spectrum");
965 }
966
967 /* Perform the normalisation of the SC vis2 and visamp
968 * to match those of the FT */
969 if (!strcmp (gravi_param_get_string (parlist, "gravity.vis.vis-correction-sc"), "FORCE")) {
970
971 cpl_msg_info (cpl_func, "Align the SC visibilities on the FT");
972 gravi_normalize_sc_to_ft (vis_data);
973
974 } else {
975 cpl_msg_info (cpl_func, "Don't align the SC visibilities on the FT");
976 }
977
978 /* Correct the wavelength due to target color shifting */
979 if (gravi_param_get_bool (parlist, "gravity.vis.color-wave-correction") ) {
980 cpl_msg_info (cpl_func,"Compute the wavelength shift due to target color");
981 gravi_wave_correct_color (vis_data);
982 CPLCHECK_CLEAN ("Cannot compute the wavelength in OI_VIS");
983 } else {
984 cpl_msg_info (cpl_func, "Don't compute the wavelength shift due to target color");
985 }
986
987 /* Co-add the observations if requested */
988 if (gravi_param_get_bool (parlist, "gravity.postprocess.average-vis")) {
989 gravi_average_vis (vis_data);
990
991 } else {
992 cpl_msg_info (cpl_func, "Don't average the different observation (if any)");
993 }
994
995 /* Recompute the TIME column from the MJD column
996 * in all OIFITS tables to follow standard */
997 gravi_vis_mjd_to_time (vis_data);
998
999 /* Save the output data file based on the first frame of the frameset in order of MJD-OBS*/
1000 const cpl_frame *it_frame;
1001 cpl_frameset * science_frames = gravi_frameset_extract_fringe_data(used_frameset);
1002 cpl_frameset_iterator *it = cpl_frameset_iterator_new(science_frames);
1003 double mjd_obs_first = DBL_MAX;
1004 while ((it_frame = cpl_frameset_iterator_get(it)) != NULL) {
1005 cpl_propertylist * this_frame_header = cpl_propertylist_load(cpl_frame_get_filename(it_frame), 0);
1006 double mjd_obs = gravi_pfits_get_mjd(this_frame_header);
1007 if (mjd_obs < mjd_obs_first)
1008 {
1009 mjd_obs_first = mjd_obs;
1010 frame = cpl_frameset_iterator_get(it);
1011 }
1012 cpl_frameset_iterator_advance(it, 1);
1013 cpl_propertylist_delete(this_frame_header);
1014 }
1015 cpl_frameset_iterator_delete(it);
1016
1017 cpl_frameset_join (used_frameset, recipe_frameset);
1018
1019 if(gravi_param_get_bool (parlist, "gravity.vis.oifits2") ) {
1020 gravi_vis_copy_fluxdata(vis_data, 1);
1021 }
1022 gravi_data_save_new (vis_data, frameset, NULL, NULL, parlist,
1023 used_frameset, frame, "gravity_vis", NULL, proCatg);
1024 cpl_frameset_delete(science_frames);
1025
1026 CPLCHECK_CLEAN ("Cannot save the VIS product");
1027
1028 /* Terminate the function */
1029 goto cleanup;
1030
1031cleanup:
1032 /* Deallocation of all variables */
1033 cpl_msg_info(cpl_func,"Memory cleanup");
1034
1035 FREELOOP (gravi_data_delete,sky_maps,nb_sky);
1036 FREE (gravi_data_delete,dark_map);
1037 FREE (gravi_data_delete,data);
1038 FREE (gravi_data_delete,preproc_data);
1039 FREE (gravi_data_delete,profile_map);
1040 FREE (gravi_data_delete,disp_map);
1041 FREE (gravi_data_delete,wave_map);
1042 FREE (gravi_data_delete,badpix_map);
1043 FREE (gravi_data_delete,p2vm_map);
1044 FREE (gravi_data_delete,p2vmred_data);
1045 FREE (gravi_data_delete,vis_data);
1046 FREE (gravi_data_delete,tmpvis_data);
1047 FREE (gravi_data_delete,diamcat_data);
1048 FREE (gravi_data_delete,static_param_data);
1049 FREE (gravi_data_delete,diodepos_data);
1050 FREE (gravi_data_delete,pca_calib_data);
1051 FREE (gravi_data_delete,eop_map);
1052 FREE (cpl_frameset_delete,darkcalib_frameset);
1053 FREE (cpl_frameset_delete,wavecalib_frameset);
1054 FREE (cpl_frameset_delete,flatcalib_frameset);
1055 FREE (cpl_frameset_delete,badcalib_frameset);
1056 FREE (cpl_frameset_delete,p2vmcalib_frameset);
1057 FREE (cpl_frameset_delete,metpos_frameset);
1058 FREE (cpl_frameset_delete,dark_frameset);
1059 FREE (cpl_frameset_delete,diamcat_frameset);
1060 FREE (cpl_frameset_delete,sky_frameset);
1061 FREE (cpl_frameset_delete,dispcalib_frameset);
1062 FREE (cpl_frameset_delete,pcacalib_frameset);
1063 FREE (cpl_frameset_delete,eop_frameset);
1064 FREE (cpl_frameset_delete,patch_frameset);
1065 FREE (cpl_frameset_delete,static_param_frameset);
1066 FREE (cpl_frameset_delete,recipe_frameset);
1067 FREE (cpl_frameset_delete,current_frameset);
1068 FREE (cpl_frameset_delete,used_frameset);
1069 FREE (cpl_free,proCatg);
1070 FREE (cpl_free,redCatg);
1071 FREE (cpl_free,skyCatg);
1072 FREE (cpl_free,mode);
1073 FREE (cpl_free,input_data_type);
1074 FREELOOP(cpl_propertylist_delete, p2vm_qcs, nb_frame);
1075
1076
1078 return (int)cpl_error_get_code();
1079}
typedefCPL_BEGIN_DECLS struct _gravi_data_ gravi_data
Definition: gravi_data.h:39
#define gravi_data_get_header(data)
Definition: gravi_data.h:75
cpl_error_code gravi_metrology_demodulate(gravi_data *met_data, cpl_boolean zero_subtracted)
Demodulate the metrology.
#define GRAVI_RECIPE_OUTPUT
Definition: gravi_dfs.h:39
#define GRAVI_P2VMRED_SINGLE_CALIB
Definition: gravi_dfs.h:63
#define GRAVI_SINGLE_SCIENCE_RAW
Definition: gravi_dfs.h:52
#define GRAVI_PREPROC
Definition: gravi_dfs.h:60
#define GRAVI_SINGLE_SKY_RAW
Definition: gravi_dfs.h:56
#define GRAVI_P2VM_MAP
Definition: gravi_dfs.h:76
#define GRAVI_P2VMRED_DUAL_CALIB
Definition: gravi_dfs.h:65
#define GRAVI_DIODE_POSITION
Definition: gravi_dfs.h:81
#define GRAVI_SINGLE_SKY_MAP
Definition: gravi_dfs.h:89
#define GRAVI_RECIPE_FLOW
Definition: gravi_dfs.h:37
#define GRAVI_VIS_DUAL_SCIENCE
Definition: gravi_dfs.h:98
#define GRAVI_VIS_SINGLE_SCIENCE
Definition: gravi_dfs.h:96
#define GRAVI_ASTROREDUCED
Definition: gravi_dfs.h:67
#define GRAVI_BAD_MAP
Definition: gravi_dfs.h:73
#define GRAVI_DUAL_CALIB_RAW
Definition: gravi_dfs.h:53
#define GRAVI_P2VMRED_SINGLE_SCIENCE
Definition: gravi_dfs.h:64
#define GRAVI_WAVE_MAP
Definition: gravi_dfs.h:75
#define GRAVI_DISP_MODEL
Definition: gravi_dfs.h:79
#define GRAVI_FLAT_MAP
Definition: gravi_dfs.h:74
#define GRAVI_DIAMETER_CAT
Definition: gravi_dfs.h:80
#define GRAVI_DARK_MAP
Definition: gravi_dfs.h:77
#define GRAVI_DUAL_SKY_MAP
Definition: gravi_dfs.h:90
#define GRAVI_RECIPE_INPUT
Definition: gravi_dfs.h:38
#define GRAVI_SINGLE_CALIB_RAW
Definition: gravi_dfs.h:51
#define GRAVI_VIS_DUAL_CALIB
Definition: gravi_dfs.h:99
#define GRAVI_DUAL_SCIENCE_RAW
Definition: gravi_dfs.h:54
#define GRAVI_SPECTRUM_ALIGNED
Definition: gravi_dfs.h:62
#define GRAVI_VIS_SINGLE_CALIB
Definition: gravi_dfs.h:97
#define GRAVI_SPECTRUM
Definition: gravi_dfs.h:61
#define GRAVI_P2VMRED_DUAL_SCIENCE
Definition: gravi_dfs.h:66
cpl_msg_info(cpl_func, "Compute WAVE_SCAN for %s", GRAVI_TYPE(type_data))
#define GRAVI_OPDC_EXT
Definition: gravi_pfits.h:62
#define GRAVI_FDDL_EXT
Definition: gravi_pfits.h:75
#define INSNAME_ACQ
Definition: gravi_pfits.h:199
@ GRAVI_DET_ALL
Definition: gravi_pfits.h:173
#define GRAVI_IMAGING_DATA_ACQ_EXT
Definition: gravi_pfits.h:41
#define GRAVI_ARRAY_GEOMETRY_EXT
Definition: gravi_pfits.h:84
#define GRAVI_METROLOGY_EXT
Definition: gravi_pfits.h:60
#define GRAVI_OPTICAL_TRAIN_EXT
Definition: gravi_pfits.h:85
#define CPLCHECK_CLEAN(msg)
Definition: gravi_utils.h:54
#define gravi_msg_function_exit(flag)
Definition: gravi_utils.h:85
#define FREE(function, variable)
Definition: gravi_utils.h:69
#define gravi_msg_function_start(flag)
Definition: gravi_utils.h:84
#define CPLCHECK_MSG(msg)
Definition: gravi_utils.h:45
#define gravi_data_check_shutter_open(data)
Definition: gravi_utils.h:91
#define FREELOOP(function, variable, n)
Definition: gravi_utils.h:72
static int gravity_vis_exec(cpl_plugin *)
Execute the plugin instance given by the interface.
Definition: gravity_vis.c:293
static char gravity_vis_description[]
Definition: gravity_vis.c:78
static int gravity_vis(cpl_frameset *, cpl_parameterlist *)
Compute the visibilities, and closure phase and create the io fits file.
Definition: gravity_vis.c:391
int cpl_plugin_get_info(cpl_pluginlist *list)
Build the list of available plugins, for this module.
Definition: gravity_vis.c:125
static char gravity_vis_short[]
Definition: gravity_vis.c:77
static int gravity_vis_create(cpl_plugin *)
Setup the recipe options
Definition: gravity_vis.c:166
static int gravity_vis_destroy(cpl_plugin *)
Destroy what has been created by the 'create' function.
Definition: gravity_vis.c:358
cpl_error_code gravi_preproc_acqcam(gravi_data *output_data, gravi_data *input_data, gravi_data *bad_map)
Preprocess the ACQ images: correct bad pixels, clean from pupil background via blinking,...
Definition: gravi_acqcam.c:208
cpl_error_code gravi_reduce_acqcam(gravi_data *output_data, gravi_data *input_data, gravi_data *sky_data, gravi_data *dark_data, gravi_data *static_param_data, int saveAcqTable)
Reduce the ACQ camera images.
cpl_error_code gravi_flatten_vis(gravi_data *vis_data, gravi_data *calib_data)
Use PCA model to flatten observed visphi. The flattened data are added to the existing VIS table.
Definition: gravi_calib.c:3279
gravi_data * gravi_average_dark(gravi_data **data, cpl_size ndata)
Average several DARK calibration map.
Definition: gravi_calib.c:457
gravi_data * gravi_compute_dark(gravi_data *raw_data)
Compute the DARK calibration map.
Definition: gravi_calib.c:125
int gravi_data_patch(gravi_data *file_to_patch, cpl_frameset *patch_frameset)
Load a RAW FITS file and create a gravi_data.
Definition: gravi_data.c:635
cpl_error_code gravi_data_clean_for_astro(gravi_data *data)
Clean the data to keep only OIFITS extensions related to SC.
Definition: gravi_data.c:2435
gravi_data * gravi_data_load_frame(cpl_frame *frame, cpl_frameset *used_frameset)
Load a FITS file and create a gravi_data.
Definition: gravi_data.c:599
cpl_error_code gravi_data_detector_cleanup(gravi_data *data, const cpl_parameterlist *parlist)
Perform self-bias correction to the SC raw data.
Definition: gravi_data.c:1232
gravi_data * gravi_data_load_rawframe(cpl_frame *frame, cpl_frameset *used_frameset)
Load a RAW FITS file and create a gravi_data.
Definition: gravi_data.c:716
cpl_error_code gravi_data_move_ext(gravi_data *output, gravi_data *input, const char *name)
Move extensions from one data to another.
Definition: gravi_data.c:1741
cpl_error_code gravi_data_append(gravi_data *first, const gravi_data *second, int force)
Append a gravi_data into another existing one.
Definition: gravi_data.c:308
cpl_error_code gravi_data_save_new(gravi_data *self, cpl_frameset *allframes, const char *filename, const char *suffix, const cpl_parameterlist *parlist, cpl_frameset *usedframes, cpl_frame *frame, const char *recipe, cpl_propertylist *applist, const char *proCatg)
Save a gravi data in a CPL-complian FITS file.
Definition: gravi_data.c:925
cpl_error_code gravi_data_copy_ext_insname(gravi_data *output, gravi_data *input, const char *name, const char *insname)
Copy extensions from one data to another.
Definition: gravi_data.c:1587
void gravi_data_delete(gravi_data *self)
Delete a gravi data.
Definition: gravi_data.c:146
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
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_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_fringe_data(cpl_frameset *frameset)
Definition: gravi_dfs.c:1320
cpl_error_code gravi_parameter_add_compute_snr(cpl_parameterlist *self)
Definition: gravi_dfs.c:776
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_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_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_parameter * gravi_parameter_add_spectrum_file(cpl_parameterlist *self)
Definition: gravi_dfs.c:506
cpl_frameset * gravi_frameset_extract_pca_calib(cpl_frameset *frameset)
Definition: gravi_dfs.c:1435
cpl_frameset * gravi_frameset_extract_dark_map(cpl_frameset *frameset)
Definition: gravi_dfs.c:1403
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_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_frameset * gravi_frameset_extract_static_param(cpl_frameset *frameset)
Definition: gravi_dfs.c:1427
cpl_error_code gravi_parameter_add_compute_vis(cpl_parameterlist *self, int isCalib)
Definition: gravi_dfs.c:940
cpl_parameter * gravi_parameter_add_biassub_file(cpl_parameterlist *self)
Definition: gravi_dfs.c:491
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
cpl_error_code gravi_compute_pointing_uv(gravi_data *p2vmred_data, gravi_data *eop_data)
Compute the pointing directions and projected baselines in OI_VIS.
Definition: gravi_eop.c:466
cpl_error_code gravi_metrology_reduce(gravi_data *data, gravi_data *eop_data, gravi_data *static_param_data, gravi_data *met_pos, const cpl_parameterlist *parlist)
Reduce the metrology.
cpl_error_code gravi_compute_qc_injection(gravi_data *data)
Compute the QC for the injection stability.
cpl_error_code gravi_compute_tau0(gravi_data *data)
Compute the QC TAU0 parameter.
cpl_error_code gravi_compute_opdc_state(gravi_data *p2vmred_data)
Compute the real-time tracking state from OPDC.
gravi_data * gravi_compute_p2vmred(gravi_data *preproc_data, gravi_data *p2vm_map, const char *mode, const cpl_parameterlist *parlist)
Converts preprocessed data into coherent fluxes using the P2VM.
cpl_error_code gravi_compute_qc_ft_opd_estimator(gravi_data *p2vmred_data)
Compute the QC for the FT linearity.
cpl_parameter * gravi_pfits_get_extrapixel_param(const cpl_propertylist *header)
Extract parameters from a product header.
Definition: gravi_pfits.c:1058
double gravi_pfits_get_mjd(const cpl_propertylist *plist)
Definition: gravi_pfits.c:526
cpl_error_code gravi_align_spectrum(gravi_data *spectrum_data, gravi_data *wave_map, gravi_data *p2vm_map, enum gravi_detector_type det_type)
Regrid the regions into a common wavelength (in-place)
gravi_data * gravi_extract_spectrum(gravi_data *raw_data, gravi_data *profile_map, gravi_data *dark_map, gravi_data *bad_map, gravi_data *sky_map, const cpl_parameterlist *parlist, enum gravi_detector_type det_type)
Create the SPECTRUM gravi_data with extracted spectrum per region.
cpl_error_code gravi_subtract_met_dark(gravi_data *preproc_data, gravi_data *dark_map)
Substract metrology dark.
cpl_error_code gravi_compute_snr(gravi_data *p2vmred_data, const cpl_parameterlist *parlist)
Compute real-time SNR and Group-Delay of the observation.
Definition: gravi_signal.c:620
cpl_error_code gravi_compute_signals(gravi_data *p2vmred_data, gravi_data *disp_data, const cpl_parameterlist *parlist)
Create intermediate signal in the P2VMREDUCED file.
cpl_error_code gravi_compute_outliers(gravi_data *p2vmred_data, const cpl_parameterlist *parlist)
Compute the outliers flags.
Definition: gravi_signal.c:536
cpl_error_code gravi_compute_rejection(gravi_data *p2vmred_data, const cpl_parameterlist *parlist)
Create rejection flags P2VMREDUCED file.
cpl_error_code gravi_copy_p2vm_qcs(gravi_data *p2vmred_data, cpl_propertylist *plist)
Copy PFACTOR and VFACTOR QCs so that they may be aggregated over all frames.
cpl_error_code gravi_compute_tf_qc(gravi_data *oi_vis, gravi_data *diamcat_data)
Fill QC parameters related to transfer function.
Definition: gravi_tf.c:1035
const char * gravi_get_license(void)
Get the pipeline copyright and license.
Definition: gravi_utils.c:104
cpl_error_code gravi_normalize_sc_to_ft(gravi_data *vis_data)
Align the SC visibilities on the FT visibilities.
Definition: gravi_vis.c:2596
gravi_data * gravi_compute_vis(gravi_data *p2vmred_data, const cpl_parameterlist *parlist, cpl_size *current_frame)
The function average the individual frames of a P2VMREDUCED file into a final, single observation per...
Definition: gravi_vis.c:1675
cpl_error_code gravi_vis_force_time(gravi_data *oi_data)
Force all data in OI_TABLE to have the same TIME and MJD.
Definition: gravi_vis.c:4163
cpl_error_code gravi_vis_mjd_to_time(gravi_data *vis_data)
Recompute the TIME column of all OIFITS extension from the MJD column, following the OIFITS standard ...
Definition: gravi_vis.c:2682
cpl_error_code gravi_vis_copy_fluxdata(gravi_data *oi_data, int delete_flux)
Duplicate the column FLUX into FLUXDATA, for OIFITS2 compliance.
Definition: gravi_vis.c:3675
cpl_error_code gravi_compute_vis_qc(gravi_data *vis_data, cpl_frameset *frameset, cpl_propertylist **frame_qcs, cpl_size nb_frame, char *input_data_type)
Compute the QC parameters for a VIS (averaged) data.
Definition: gravi_vis.c:2237
cpl_error_code gravi_flat_flux(gravi_data *vis_data, gravi_data *p2vm_map)
Divide the OI_FLUX by OI_FLUX from the P2VM (no checks, no time distance...)
Definition: gravi_vis.c:2736
cpl_error_code gravi_average_vis(gravi_data *oi_data)
Coadd the observations together.
Definition: gravi_vis.c:3030
cpl_error_code gravi_wave_correct_color(gravi_data *vis_data)
Create a OI_WAVELENGTH_CORR table with color corrected wavelength.
Definition: gravi_wave.c:1748