FORS Pipeline Reference Manual 4.9.9
fors_science.c
00001 /* $Id: fors_science.c,v 1.31 2011/12/01 15:32:26 cgarcia Exp $
00002  *
00003  * This file is part of the FORS Data Reduction Pipeline
00004  * Copyright (C) 2002-2010 European Southern Observatory
00005  *
00006  * This program is free software; you can redistribute it and/or modify
00007  * it under the terms of the GNU General Public License as published by
00008  * the Free Software Foundation; either version 2 of the License, or
00009  * (at your option) any later version.
00010  *
00011  * This program is distributed in the hope that it will be useful,
00012  * but WITHOUT ANY WARRANTY; without even the implied warranty of
00013  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00014  * GNU General Public License for more details.
00015  *
00016  * You should have received a copy of the GNU General Public License
00017  * along with this program; if not, write to the Free Software
00018  * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
00019  */
00020 
00021 /*
00022  * $Author: cgarcia $
00023  * $Date: 2011/12/01 15:32:26 $
00024  * $Revision: 1.31 $
00025  * $Name: fors-4_9_9 $
00026  */
00027 
00028 #ifdef HAVE_CONFIG_H
00029 #include <config.h>
00030 #endif
00031 
00032 #include <math.h>
00033 #include <string.h>
00034 #include <cpl.h>
00035 #include <moses.h>
00036 #include <fors_dfs.h>
00037 #include <fors_tools.h>
00038 #include <fors_qc.h>
00039 #include <fors_utils.h>
00040 
00041 static int fors_science_create(cpl_plugin *);
00042 static int fors_science_exec(cpl_plugin *);
00043 static int fors_science_destroy(cpl_plugin *);
00044 static int fors_science(cpl_parameterlist *, cpl_frameset *);
00045 
00046 static char fors_science_description[] =
00047 "This recipe is used to reduce scientific spectra using the extraction\n"
00048 "mask and the products created by the recipe fors_calib. The spectra are\n"
00049 "bias subtracted, flat fielded (if a normalised flat field is specified)\n"
00050 "and remapped eliminating the optical distortions. The wavelength calibration\n"
00051 "can be optionally upgraded using a number of sky lines: if no sky lines\n"
00052 "catalog of wavelengths is specified, an internal one is used instead.\n"
00053 "If the alignment to the sky lines is performed, the input dispersion\n"
00054 "coefficients table is upgraded and saved to disk, and a new CCD wavelengths\n"
00055 "map is created.\n"
00056 "This recipe accepts both FORS1 and FORS2 frames. A grism table (typically\n"
00057 "depending on the instrument mode, and in particular on the grism used)\n"
00058 "may also be specified: this table contains a default recipe parameter\n" 
00059 "setting to control the way spectra are extracted for a specific instrument\n"
00060 "mode, as it is used for automatic run of the pipeline on Paranal and in\n" 
00061 "Garching. If this table is specified, it will modify the default recipe\n" 
00062 "parameter setting, with the exception of those parameters which have been\n" 
00063 "explicitly modifyed on the command line. If a grism table is not specified,\n"
00064 "the input recipe parameters values will always be read from the command\n" 
00065 "line, or from an esorex configuration file if present, or from their\n" 
00066 "generic default values (that are rarely meaningful).\n" 
00067 "In the table below the MXU acronym can be read alternatively as MOS\n"
00068 "and LSS, depending on the instrument mode of the input data. The acronym\n"
00069 "SCI on products should be read STD in case of standard stars observations\n"
00070 "A CURV_COEFF table is not (yet) expected for LSS data.\n"
00071 "Either a scientific or a standard star exposure can be specified in input.\n"
00072 "Only in case of a standard star exposure input, the atmospheric extinction\n"
00073 "table and a table with the physical fluxes of the observed standard star\n"
00074 "must be specified in input, and a spectro-photometric table is created in\n"
00075 "output. This table can then be input again to this recipe, always with an\n"
00076 "atmospheric extinction table, and if a photometric calibration is requested\n"
00077 "then flux calibrated spectra (in units of erg/cm/cm/s/Angstrom) are also\n" 
00078 "written in output.\n\n"
00079 
00080 "Input files:\n\n"
00081 "  DO category:               Type:       Explanation:         Required:\n"
00082 "  SCIENCE_MXU                Raw         Scientific exposure     Y\n"
00083 "  or STANDARD_MXU            Raw         Standard star exposure  Y\n"
00084 "  MASTER_BIAS                Calib       Master bias             Y\n"
00085 "  GRISM_TABLE                Calib       Grism table             .\n"
00086 "  MASTER_SKYLINECAT          Calib       Sky lines catalog       .\n"
00087 "\n"
00088 "  MASTER_NORM_FLAT_MXU       Calib       Normalised flat field   .\n"
00089 "  DISP_COEFF_MXU             Calib       Inverse dispersion      Y\n"
00090 "  CURV_COEFF_MXU             Calib       Spectral curvature      Y\n"
00091 "  SLIT_LOCATION_MXU          Calib       Slits positions table   Y\n"
00092 "\n"
00093 "  or, in case of LSS-like MOS/MXU data,\n"
00094 "\n"
00095 "  MASTER_NORM_FLAT_LONG_MXU  Calib       Normalised flat field   .\n"
00096 "  DISP_COEFF_LONG_MXU        Calib       Inverse dispersion      Y\n"
00097 "\n"
00098 "  In case STANDARD_MXU is specified in input,\n"
00099 "\n"
00100 "  EXTINCT_TABLE              Calib       Atmospheric extinction  Y\n"
00101 "  STD_FLUX_TABLE             Calib       Standard star flux      Y\n"
00102 "\n"
00103 "  In case a photometric calibration is requested for scientific\n"
00104 "  data, the following inputs are mandatory:\n"
00105 "\n"
00106 "  EXTINCT_TABLE              Calib       Atmospheric extinction  Y\n"
00107 "  SPECPHOT_TABLE             Calib       Response curves         Y\n"
00108 "\n"
00109 "  If requested for standard star data, the SPECPHOT_TABLE can be dropped:\n"
00110 "  in this case the correction is applied using the SPECPHOT_TABLE produced\n"
00111 "  in the same run.\n"
00112 "\n"
00113 "Output files:\n"
00114 "\n"
00115 "  DO category:               Data type:  Explanation:\n"
00116 "  REDUCED_SCI_MXU            FITS image  Extracted scientific spectra\n"
00117 "  REDUCED_SKY_SCI_MXU        FITS image  Extracted sky spectra\n"
00118 "  REDUCED_ERROR_SCI_MXU      FITS image  Errors on extracted spectra\n"
00119 "  UNMAPPED_SCI_MXU           FITS image  Sky subtracted scientific spectra\n"
00120 "  MAPPED_SCI_MXU             FITS image  Rectified scientific spectra\n"
00121 "  MAPPED_ALL_SCI_MXU         FITS image  Rectified science spectra with sky\n"
00122 "  MAPPED_SKY_SCI_MXU         FITS image  Rectified sky spectra\n"
00123 "  UNMAPPED_SKY_SCI_MXU       FITS image  Sky on CCD\n"
00124 "  GLOBAL_SKY_SPECTRUM_MXU    FITS table  Global sky spectrum\n"
00125 "  OBJECT_TABLE_SCI_MXU       FITS table  Positions of detected objects\n"
00126 "\n"
00127 "  Only if the sky-alignment of the wavelength solution is requested:\n"
00128 "  SKY_SHIFTS_LONG_SCI_MXU    FITS table  Sky lines offsets (LSS-like data)\n"
00129 "  or SKY_SHIFTS_SLIT_SCI_MXU FITS table  Sky lines offsets (MOS-like data)\n"
00130 "  DISP_COEFF_SCI_MXU         FITS table  Upgraded dispersion coefficients\n"
00131 "  WAVELENGTH_MAP_SCI_MXU     FITS image  Upgraded wavelength map\n"
00132 "\n"
00133 "  Only if a STANDARD_MXU is specified in input:\n"
00134 "  SPECPHOT_TABLE             FITS table  Efficiency and response curves\n"
00135 "\n"
00136 "  Only if a photometric calibration was requested:\n"
00137 "  REDUCED_FLUX_SCI_MXU       FITS image  Flux calibrated scientific spectra\n"
00138 "  REDUCED_FLUX_ERROR_SCI_MXU FITS image  Errors on flux calibrated spectra\n"
00139 "  MAPPED_FLUX_SCI_MXU        FITS image  Flux calibrated slit spectra\n\n";
00140 
00141 #define fors_science_exit(message)            \
00142 {                                             \
00143 if (message) cpl_msg_error(recipe, message);  \
00144 cpl_free(exptime);                            \
00145 cpl_free(instrume);                           \
00146 cpl_image_delete(dummy);                      \
00147 cpl_image_delete(mapped);                     \
00148 cpl_image_delete(mapped_sky);                 \
00149 cpl_image_delete(mapped_cleaned);             \
00150 cpl_image_delete(skylocalmap);                \
00151 cpl_image_delete(skymap);                     \
00152 cpl_image_delete(smapped);                    \
00153 cpl_table_delete(offsets);                    \
00154 cpl_table_delete(photcal);                    \
00155 cpl_table_delete(sky);                        \
00156 cpl_image_delete(bias);                       \
00157 cpl_image_delete(spectra);                    \
00158 cpl_image_delete(coordinate);                 \
00159 cpl_image_delete(norm_flat);                  \
00160 cpl_image_delete(rainbow);                    \
00161 cpl_image_delete(rectified);                  \
00162 cpl_image_delete(wavemap);                    \
00163 cpl_image_delete(wavemaplss);                 \
00164 cpl_propertylist_delete(qclist);              \
00165 cpl_propertylist_delete(header);              \
00166 cpl_propertylist_delete(save_header);         \
00167 cpl_table_delete(grism_table);                \
00168 cpl_table_delete(idscoeff);                   \
00169 cpl_table_delete(maskslits);                  \
00170 cpl_table_delete(overscans);                  \
00171 cpl_table_delete(polytraces);                 \
00172 cpl_table_delete(slits);                      \
00173 cpl_table_delete(wavelengths);                \
00174 cpl_vector_delete(lines);                     \
00175 cpl_msg_indent_less();                        \
00176 return -1;                                    \
00177 }
00178 
00179 
00180 #define fors_science_exit_memcheck(message)   \
00181 {                                             \
00182 if (message) cpl_msg_info(recipe, message);   \
00183 cpl_free(exptime);                            \
00184 cpl_free(instrume);                           \
00185 cpl_image_delete(dummy);                      \
00186 cpl_image_delete(mapped);                     \
00187 cpl_image_delete(mapped_cleaned);             \
00188 cpl_image_delete(mapped_sky);                 \
00189 cpl_image_delete(skylocalmap);                \
00190 cpl_image_delete(skymap);                     \
00191 cpl_image_delete(smapped);                    \
00192 cpl_table_delete(offsets);                    \
00193 cpl_table_delete(photcal);                    \
00194 cpl_table_delete(sky);                        \
00195 cpl_image_delete(bias);                       \
00196 cpl_image_delete(spectra);                    \
00197 cpl_image_delete(coordinate);                 \
00198 cpl_image_delete(norm_flat);                  \
00199 cpl_image_delete(rainbow);                    \
00200 cpl_image_delete(rectified);                  \
00201 cpl_image_delete(wavemap);                    \
00202 cpl_image_delete(wavemaplss);                 \
00203 cpl_propertylist_delete(header);              \
00204 cpl_propertylist_delete(save_header);         \
00205 cpl_table_delete(grism_table);                \
00206 cpl_table_delete(idscoeff);                   \
00207 cpl_table_delete(maskslits);                  \
00208 cpl_table_delete(overscans);                  \
00209 cpl_table_delete(polytraces);                 \
00210 cpl_table_delete(slits);                      \
00211 cpl_table_delete(wavelengths);                \
00212 cpl_vector_delete(lines);                     \
00213 cpl_msg_indent_less();                        \
00214 return 0;                                     \
00215 }
00216 
00217 
00229 int cpl_plugin_get_info(cpl_pluginlist *list)
00230 {
00231     cpl_recipe *recipe = cpl_calloc(1, sizeof *recipe );
00232     cpl_plugin *plugin = &recipe->interface;
00233 
00234     cpl_plugin_init(plugin,
00235                     CPL_PLUGIN_API,
00236                     FORS_BINARY_VERSION,
00237                     CPL_PLUGIN_TYPE_RECIPE,
00238                     "fors_science",
00239                     "Extraction of scientific spectra",
00240                     fors_science_description,
00241                     "Carlo Izzo",
00242                     PACKAGE_BUGREPORT,
00243                     fors_get_license(),
00244                     fors_science_create,
00245                     fors_science_exec,
00246                     fors_science_destroy);
00247 
00248     cpl_pluginlist_append(list, plugin);
00249     
00250     return 0;
00251 }
00252 
00253 
00264 static int fors_science_create(cpl_plugin *plugin)
00265 {
00266     cpl_recipe    *recipe;
00267     cpl_parameter *p;
00268 
00269 
00270     /* 
00271      * Check that the plugin is part of a valid recipe 
00272      */
00273 
00274     if (cpl_plugin_get_type(plugin) == CPL_PLUGIN_TYPE_RECIPE) 
00275         recipe = (cpl_recipe *)plugin;
00276     else 
00277         return -1;
00278 
00279     /* 
00280      * Create the parameters list in the cpl_recipe object 
00281      */
00282 
00283     recipe->parameters = cpl_parameterlist_new(); 
00284 
00285 
00286     /*
00287      * Dispersion
00288      */
00289 
00290     p = cpl_parameter_new_value("fors.fors_science.dispersion",
00291                                 CPL_TYPE_DOUBLE,
00292                                 "Resampling step (Angstrom/pixel)",
00293                                 "fors.fors_science",
00294                                 0.0);
00295     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "dispersion");
00296     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00297     cpl_parameterlist_append(recipe->parameters, p);
00298 
00299     /*
00300      * Sky lines alignment
00301      */
00302 
00303     p = cpl_parameter_new_value("fors.fors_science.skyalign",
00304                                 CPL_TYPE_INT,
00305                                 "Polynomial order for sky lines alignment, "
00306                                 "or -1 to avoid alignment",
00307                                 "fors.fors_science",
00308                                 0);
00309     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "skyalign");
00310     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00311     cpl_parameterlist_append(recipe->parameters, p);
00312 
00313     /*
00314      * Line catalog table column containing the sky reference wavelengths
00315      */
00316 
00317     p = cpl_parameter_new_value("fors.fors_science.wcolumn",
00318                                 CPL_TYPE_STRING,
00319                                 "Name of sky line catalog table column "
00320                                 "with wavelengths",
00321                                 "fors.fors_science",
00322                                 "WLEN");
00323     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "wcolumn");
00324     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00325     cpl_parameterlist_append(recipe->parameters, p);
00326 
00327     /*
00328      * Start wavelength for spectral extraction
00329      */
00330 
00331     p = cpl_parameter_new_value("fors.fors_science.startwavelength",
00332                                 CPL_TYPE_DOUBLE,
00333                                 "Start wavelength in spectral extraction",
00334                                 "fors.fors_science",
00335                                 0.0);
00336     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "startwavelength");
00337     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00338     cpl_parameterlist_append(recipe->parameters, p);
00339 
00340     /*
00341      * End wavelength for spectral extraction
00342      */
00343 
00344     p = cpl_parameter_new_value("fors.fors_science.endwavelength",
00345                                 CPL_TYPE_DOUBLE,
00346                                 "End wavelength in spectral extraction",
00347                                 "fors.fors_science",
00348                                 0.0);
00349     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "endwavelength");
00350     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00351     cpl_parameterlist_append(recipe->parameters, p);
00352 
00353     /*
00354      * Flux conservation
00355      */
00356 
00357     p = cpl_parameter_new_value("fors.fors_science.flux",
00358                                 CPL_TYPE_BOOL,
00359                                 "Apply flux conservation",
00360                                 "fors.fors_science",
00361                                 TRUE);
00362     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "flux");
00363     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00364     cpl_parameterlist_append(recipe->parameters, p);
00365 
00366     /*
00367      * Apply flat field
00368      */
00369 
00370     p = cpl_parameter_new_value("fors.fors_science.flatfield",
00371                                 CPL_TYPE_BOOL,
00372                                 "Apply flat field",
00373                                 "fors.fors_science",
00374                                 TRUE);
00375     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "flatfield");
00376     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00377     cpl_parameterlist_append(recipe->parameters, p);
00378 
00379     /*
00380      * Global sky subtraction
00381      */
00382 
00383     p = cpl_parameter_new_value("fors.fors_science.skyglobal",
00384                                 CPL_TYPE_BOOL,
00385                                 "Subtract global sky spectrum from CCD",
00386                                 "fors.fors_science",
00387                                 FALSE);
00388     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "skyglobal");
00389     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00390     cpl_parameterlist_append(recipe->parameters, p);
00391 
00392     /*
00393      * Local sky subtraction on extracted spectra
00394      */
00395 
00396 /*** New sky subtraction (search NSS)
00397     p = cpl_parameter_new_value("fors.fors_science.skymedian",
00398                                 CPL_TYPE_INT,
00399                                 "Degree of sky fitting polynomial for "
00400                                 "sky subtraction from extracted "
00401                                 "slit spectra (MOS/MXU only, -1 to disable it)",
00402                                 "fors.fors_science",
00403                                 0);
00404     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "skymedian");
00405     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00406     cpl_parameterlist_append(recipe->parameters, p);
00407 ***/
00408 
00409     p = cpl_parameter_new_value("fors.fors_science.skymedian",
00410                                 CPL_TYPE_BOOL,
00411                                 "Sky subtraction from extracted slit spectra",
00412                                 "fors.fors_science",
00413                                 FALSE);
00414     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "skymedian");
00415     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00416     cpl_parameterlist_append(recipe->parameters, p);
00417 
00418     /*
00419      * Local sky subtraction on CCD spectra
00420      */
00421 
00422     p = cpl_parameter_new_value("fors.fors_science.skylocal",
00423                                 CPL_TYPE_BOOL,
00424                                 "Sky subtraction from CCD slit spectra",
00425                                 "fors.fors_science",
00426                                 TRUE);
00427     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "skylocal");
00428     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00429     cpl_parameterlist_append(recipe->parameters, p);
00430 
00431     /*
00432      * Cosmic rays removal
00433      */
00434 
00435     p = cpl_parameter_new_value("fors.fors_science.cosmics",
00436                                 CPL_TYPE_BOOL,
00437                                 "Eliminate cosmic rays hits (only if global "
00438                                 "sky subtraction is also requested)",
00439                                 "fors.fors_science",
00440                                 FALSE);
00441     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "cosmics");
00442     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00443     cpl_parameterlist_append(recipe->parameters, p);
00444 
00445     /*
00446      * Slit margin
00447      */
00448 
00449     p = cpl_parameter_new_value("fors.fors_science.slit_margin",
00450                                 CPL_TYPE_INT,
00451                                 "Number of pixels to exclude at each slit "
00452                                 "in object detection and extraction",
00453                                 "fors.fors_science",
00454                                 3);
00455     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "slit_margin");
00456     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00457     cpl_parameterlist_append(recipe->parameters, p);
00458 
00459     /*
00460      * Extraction radius
00461      */
00462 
00463     p = cpl_parameter_new_value("fors.fors_science.ext_radius",
00464                                 CPL_TYPE_INT,
00465                                 "Maximum extraction radius for detected "
00466                                 "objects (pixel)",
00467                                 "fors.fors_science",
00468                                 12);
00469     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "ext_radius");
00470     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00471     cpl_parameterlist_append(recipe->parameters, p);
00472 
00473     /*
00474      * Contamination radius
00475      */
00476 
00477     p = cpl_parameter_new_value("fors.fors_science.cont_radius",
00478                                 CPL_TYPE_INT,
00479                                 "Minimum distance at which two objects "
00480                                 "of equal luminosity do not contaminate "
00481                                 "each other (pixel)",
00482                                 "fors.fors_science",
00483                                 0);
00484     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "cont_radius");
00485     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00486     cpl_parameterlist_append(recipe->parameters, p);
00487 
00488     /*
00489      * Object extraction method
00490      */
00491 
00492     p = cpl_parameter_new_value("fors.fors_science.ext_mode",
00493                                 CPL_TYPE_INT,
00494                                 "Object extraction method: 0 = aperture, "
00495                                 "1 = Horne optimal extraction",
00496                                 "fors.fors_science",
00497                                 1);
00498     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "ext_mode");
00499     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00500     cpl_parameterlist_append(recipe->parameters, p);
00501 
00502     /*
00503      * Normalise output by exposure time
00504      */
00505 
00506     p = cpl_parameter_new_value("fors.fors_science.time_normalise",
00507                                 CPL_TYPE_BOOL,
00508                                 "Normalise output spectra by the exposure time",
00509                                 "fors.fors_science",
00510                                 TRUE);
00511     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "time_normalise");
00512     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00513     cpl_parameterlist_append(recipe->parameters, p);
00514 
00515     /*
00516      * Order of polynomial modeling the instrument response.
00517      */
00518 
00519     p = cpl_parameter_new_value("fors.fors_science.response",
00520                                 CPL_TYPE_INT,
00521                                 "Order of polynomial modeling the "
00522                                 "instrument response",
00523                                 "fors.fors_science",
00524                                 5);
00525     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "response");
00526     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00527     cpl_parameterlist_append(recipe->parameters, p);
00528 
00529     /*
00530      * Order of polynomial modeling the instrument response.
00531      */
00532 
00533     p = cpl_parameter_new_value("fors.fors_science.photometry",
00534                                 CPL_TYPE_BOOL,
00535                                 "Apply spectrophotometric calibration",
00536                                 "fors.fors_science",
00537                                 FALSE);
00538     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "photometry");
00539     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00540     cpl_parameterlist_append(recipe->parameters, p);
00541 
00542     /*
00543      * Computation of QC1 parameters
00544      */
00545 
00546     p = cpl_parameter_new_value("fors.fors_science.qc",
00547                                 CPL_TYPE_BOOL,
00548                                 "Compute QC1 parameters",
00549                                 "fors.fors_science",
00550                                 TRUE);
00551     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "qc");
00552     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00553     cpl_parameterlist_append(recipe->parameters, p);
00554 
00555     return 0;
00556 }
00557 
00558 
00567 static int fors_science_exec(cpl_plugin *plugin)
00568 {
00569     cpl_recipe *recipe;
00570     
00571     if (cpl_plugin_get_type(plugin) == CPL_PLUGIN_TYPE_RECIPE) 
00572         recipe = (cpl_recipe *)plugin;
00573     else 
00574         return -1;
00575 
00576     return fors_science(recipe->parameters, recipe->frames);
00577 }
00578 
00579 
00588 static int fors_science_destroy(cpl_plugin *plugin)
00589 {
00590     cpl_recipe *recipe;
00591     
00592     if (cpl_plugin_get_type(plugin) == CPL_PLUGIN_TYPE_RECIPE) 
00593         recipe = (cpl_recipe *)plugin;
00594     else 
00595         return -1;
00596 
00597     cpl_parameterlist_delete(recipe->parameters); 
00598 
00599     return 0;
00600 }
00601 
00602 
00612 static int fors_science(cpl_parameterlist *parlist, cpl_frameset *frameset)
00613 {
00614 
00615     const char *recipe = "fors_science";
00616 
00617 
00618     /*
00619      * Input parameters
00620      */
00621 
00622     double      dispersion;
00623     int         skyalign;
00624     const char *wcolumn;
00625     double      startwavelength;
00626     double      endwavelength;
00627     int         flux;
00628     int         flatfield;
00629     int         skyglobal;
00630     int         skylocal;
00631     int         skymedian;
00632     int         cosmics;
00633     int         slit_margin;
00634     int         ext_radius;
00635     int         cont_radius;
00636     int         ext_mode;
00637     int         time_normalise;
00638     int         res_order;
00639     int         photometry;
00640     int         qc;
00641 
00642 
00643     /*
00644      * CPL objects
00645      */
00646 
00647     cpl_imagelist    *all_science;
00648     cpl_image       **images;
00649 
00650     cpl_image        *bias           = NULL;
00651     cpl_image        *norm_flat      = NULL;
00652     cpl_image        *spectra        = NULL;
00653     cpl_image        *rectified      = NULL;
00654     cpl_image        *coordinate     = NULL;
00655     cpl_image        *rainbow        = NULL;
00656     cpl_image        *mapped         = NULL;
00657     cpl_image        *mapped_sky     = NULL;
00658     cpl_image        *mapped_cleaned = NULL;
00659     cpl_image        *smapped        = NULL;
00660     cpl_image        *wavemap        = NULL;
00661     cpl_image        *wavemaplss     = NULL;
00662     cpl_image        *skymap         = NULL;
00663     cpl_image        *skylocalmap    = NULL;
00664     cpl_image        *dummy          = NULL;
00665 
00666     cpl_table        *grism_table    = NULL;
00667     cpl_table        *overscans      = NULL;
00668     cpl_table        *wavelengths    = NULL;
00669     cpl_table        *idscoeff       = NULL;
00670     cpl_table        *slits          = NULL;
00671     cpl_table        *maskslits      = NULL;
00672     cpl_table        *polytraces     = NULL;
00673     cpl_table        *offsets        = NULL;
00674     cpl_table        *sky            = NULL;
00675     cpl_table        *photcal        = NULL;
00676 
00677     cpl_vector       *lines          = NULL;
00678 
00679     cpl_propertylist *header         = NULL;
00680     cpl_propertylist *save_header    = NULL;
00681     cpl_propertylist *qclist         = NULL;
00682 
00683 
00684     /*
00685      * Auxiliary variables
00686      */
00687 
00688     char        version[80];
00689     char       *instrume = NULL;
00690     char       *wheel4;
00691     const char *science_tag;
00692     const char *master_norm_flat_tag;
00693     const char *disp_coeff_tag;
00694     const char *disp_coeff_sky_tag;
00695     const char *wavelength_map_sky_tag;
00696     const char *curv_coeff_tag;
00697     const char *slit_location_tag;
00698     const char *reduced_science_tag;
00699     const char *reduced_flux_science_tag;
00700     const char *reduced_sky_tag;
00701     const char *reduced_error_tag;
00702     const char *reduced_flux_error_tag;
00703     const char *mapped_science_tag;
00704     const char *mapped_flux_science_tag;
00705     const char *unmapped_science_tag;
00706     const char *mapped_science_sky_tag;
00707     const char *mapped_sky_tag;
00708     const char *unmapped_sky_tag;
00709     const char *global_sky_spectrum_tag;
00710     const char *object_table_tag;
00711     const char *skylines_offsets_tag;
00712     const char *specphot_tag;
00713     const char *master_specphot_tag = "MASTER_SPECPHOT_TABLE";
00714     int         mxu, mos, lss;
00715     int         treat_as_lss = 0;
00716     int         have_phot = 0;
00717     int         nslits;
00718     int         nscience;
00719     double     *xpos;
00720     double     *exptime = NULL;
00721     double      alltime;
00722     double      airmass;
00723     double      mxpos;
00724     double      mean_rms;
00725     int         nlines;
00726     int         rebin;
00727     double     *line;
00728     int         nx, ny;
00729     int         ccd_xsize, ccd_ysize;
00730     double      reference;
00731     double      gain;
00732     double      ron;
00733     int         standard;
00734     int         highres;
00735     int         i;
00736     double      wstart;
00737     double      wstep;
00738     int         wcount;
00739 
00740 
00741     snprintf(version, 80, "%s-%s", PACKAGE, PACKAGE_VERSION);
00742 
00743     cpl_msg_set_indentation(2);
00744 
00745     if (dfs_files_dont_exist(frameset))
00746         fors_science_exit(NULL);
00747 
00748 
00749     /* 
00750      * Get configuration parameters
00751      */
00752 
00753     cpl_msg_info(recipe, "Recipe %s configuration parameters:", recipe);
00754     cpl_msg_indent_more();
00755 
00756     if (cpl_frameset_count_tags(frameset, "GRISM_TABLE") > 1)
00757         fors_science_exit("Too many in input: GRISM_TABLE");
00758 
00759     grism_table = dfs_load_table(frameset, "GRISM_TABLE", 1);
00760 
00761     dispersion = dfs_get_parameter_double(parlist, 
00762                     "fors.fors_science.dispersion", grism_table);
00763 
00764     if (dispersion <= 0.0)
00765         fors_science_exit("Invalid resampling step");
00766 
00767     skyalign = dfs_get_parameter_int(parlist, 
00768                     "fors.fors_science.skyalign", NULL);
00769 
00770     if (skyalign > 2)
00771         fors_science_exit("Max polynomial degree for sky alignment is 2");
00772 
00773     wcolumn = dfs_get_parameter_string(parlist, 
00774                     "fors.fors_science.wcolumn", NULL);
00775 
00776     startwavelength = dfs_get_parameter_double(parlist, 
00777                     "fors.fors_science.startwavelength", grism_table);
00778     if (startwavelength < 3000.0 || startwavelength > 13000.0)
00779         fors_science_exit("Invalid wavelength");
00780 
00781     endwavelength = dfs_get_parameter_double(parlist, 
00782                     "fors.fors_science.endwavelength", grism_table);
00783     if (endwavelength < 3000.0 || endwavelength > 13000.0)
00784         fors_science_exit("Invalid wavelength");
00785 
00786     if (endwavelength - startwavelength <= 0.0)
00787         fors_science_exit("Invalid wavelength interval");
00788 
00789     flux = dfs_get_parameter_bool(parlist, "fors.fors_science.flux", NULL);
00790 
00791     flatfield = dfs_get_parameter_bool(parlist, "fors.fors_science.flatfield", 
00792                                        NULL);
00793 
00794     skyglobal = dfs_get_parameter_bool(parlist, "fors.fors_science.skyglobal", 
00795                                        NULL);
00796     skylocal  = dfs_get_parameter_bool(parlist, "fors.fors_science.skylocal", 
00797                                        NULL);
00798     skymedian = dfs_get_parameter_bool(parlist, "fors.fors_science.skymedian", 
00799                                        NULL);
00800 /* NSS
00801     skymedian = dfs_get_parameter_int(parlist, "fors.fors_science.skymedian", 
00802                                        NULL);
00803 */
00804 
00805     if (skylocal && skyglobal)
00806         fors_science_exit("Cannot apply both local and global sky subtraction");
00807 
00808     if (skylocal && skymedian)
00809         fors_science_exit("Cannot apply sky subtraction both on extracted "
00810                           "and non-extracted spectra");
00811 
00812     cosmics = dfs_get_parameter_bool(parlist, 
00813                                      "fors.fors_science.cosmics", NULL);
00814 
00815     if (cosmics)
00816         if (!(skyglobal || skylocal))
00817             fors_science_exit("Cosmic rays correction requires "
00818                               "either skylocal=true or skyglobal=true");
00819 
00820     slit_margin = dfs_get_parameter_int(parlist, 
00821                                         "fors.fors_science.slit_margin",
00822                                         NULL);
00823     if (slit_margin < 0)
00824         fors_science_exit("Value must be zero or positive");
00825 
00826     ext_radius = dfs_get_parameter_int(parlist, "fors.fors_science.ext_radius",
00827                                        NULL);
00828     if (ext_radius < 0)
00829         fors_science_exit("Value must be zero or positive");
00830 
00831     cont_radius = dfs_get_parameter_int(parlist, 
00832                                         "fors.fors_science.cont_radius",
00833                                         NULL);
00834     if (cont_radius < 0)
00835         fors_science_exit("Value must be zero or positive");
00836 
00837     ext_mode = dfs_get_parameter_int(parlist, "fors.fors_science.ext_mode",
00838                                        NULL);
00839     if (ext_mode < 0 || ext_mode > 1)
00840         fors_science_exit("Invalid object extraction mode");
00841 
00842     time_normalise = dfs_get_parameter_bool(parlist, 
00843                              "fors.fors_science.time_normalise", NULL);
00844 
00845     res_order = dfs_get_parameter_int(parlist, "fors.fors_science.response",
00846                                       NULL);
00847     if (res_order < 2 || res_order > 10)
00848         fors_science_exit("Invalid instrument response modeling polynomial");
00849 
00850     photometry = dfs_get_parameter_bool(parlist,
00851                              "fors.fors_science.photometry", NULL);
00852 
00853     qc = dfs_get_parameter_bool(parlist, "fors.fors_science.qc", NULL);
00854 
00855     cpl_table_delete(grism_table); grism_table = NULL;
00856 
00857     if (cpl_error_get_code())
00858         fors_science_exit("Failure getting the configuration parameters");
00859 
00860     
00861     /* 
00862      * Check input set-of-frames
00863      */
00864 
00865     cpl_msg_indent_less();
00866     cpl_msg_info(recipe, "Check input set-of-frames:");
00867     cpl_msg_indent_more();
00868 
00869     if (!dfs_equal_keyword(frameset, "ESO INS GRIS1 ID"))
00870         fors_science_exit("Input frames are not from the same grism");
00871 
00872     if (!dfs_equal_keyword(frameset, "ESO INS FILT1 ID"))
00873         fors_science_exit("Input frames are not from the same filter");
00874 
00875     if (!dfs_equal_keyword(frameset, "ESO DET CHIP1 ID"))
00876         fors_science_exit("Input frames are not from the same chip");
00877 
00878     mxu = cpl_frameset_count_tags(frameset, "SCIENCE_MXU");
00879     mos = cpl_frameset_count_tags(frameset, "SCIENCE_MOS");
00880     lss = cpl_frameset_count_tags(frameset, "SCIENCE_LSS");
00881     standard = 0;
00882 
00883     if (mxu + mos + lss == 0) {
00884         mxu = cpl_frameset_count_tags(frameset, "STANDARD_MXU");
00885         mos = cpl_frameset_count_tags(frameset, "STANDARD_MOS");
00886         lss = cpl_frameset_count_tags(frameset, "STANDARD_LSS");
00887         standard = 1;
00888     }
00889 
00890     if (mxu + mos + lss == 0)
00891         fors_science_exit("Missing input scientific frame");
00892 
00893     nscience = mxu + mos + lss;
00894 
00895     if (mxu && mxu < nscience)
00896         fors_science_exit("Input scientific frames must be of the same type"); 
00897 
00898     if (mos && mos < nscience)
00899         fors_science_exit("Input scientific frames must be of the same type"); 
00900 
00901     if (lss && lss < nscience)
00902         fors_science_exit("Input scientific frames must be of the same type"); 
00903 
00904     if (mxu) {
00905         cpl_msg_info(recipe, "MXU data found");
00906         if (standard) {
00907             science_tag              = "STANDARD_MXU";
00908             reduced_science_tag      = "REDUCED_STD_MXU";
00909             reduced_flux_science_tag = "REDUCED_FLUX_STD_MXU";
00910             unmapped_science_tag     = "UNMAPPED_STD_MXU";
00911             mapped_science_tag       = "MAPPED_STD_MXU";
00912             mapped_flux_science_tag  = "MAPPED_FLUX_STD_MXU";
00913             mapped_science_sky_tag   = "MAPPED_ALL_STD_MXU";
00914             skylines_offsets_tag     = "SKY_SHIFTS_SLIT_STD_MXU";
00915             wavelength_map_sky_tag   = "WAVELENGTH_MAP_STD_MXU";
00916             disp_coeff_sky_tag       = "DISP_COEFF_STD_MXU";
00917             mapped_sky_tag           = "MAPPED_SKY_STD_MXU";
00918             unmapped_sky_tag         = "UNMAPPED_SKY_STD_MXU";
00919             object_table_tag         = "OBJECT_TABLE_STD_MXU";
00920             reduced_sky_tag          = "REDUCED_SKY_STD_MXU";
00921             reduced_error_tag        = "REDUCED_ERROR_STD_MXU";
00922             reduced_flux_error_tag   = "REDUCED_FLUX_ERROR_STD_MXU";
00923             specphot_tag             = "SPECPHOT_TABLE";
00924         }
00925         else {
00926             science_tag              = "SCIENCE_MXU";
00927             reduced_science_tag      = "REDUCED_SCI_MXU";
00928             reduced_flux_science_tag = "REDUCED_FLUX_SCI_MXU";
00929             unmapped_science_tag     = "UNMAPPED_SCI_MXU";
00930             mapped_science_tag       = "MAPPED_SCI_MXU";
00931             mapped_flux_science_tag  = "MAPPED_FLUX_SCI_MXU";
00932             mapped_science_sky_tag   = "MAPPED_ALL_SCI_MXU";
00933             skylines_offsets_tag     = "SKY_SHIFTS_SLIT_SCI_MXU";
00934             wavelength_map_sky_tag   = "WAVELENGTH_MAP_SCI_MXU";
00935             disp_coeff_sky_tag       = "DISP_COEFF_SCI_MXU";
00936             mapped_sky_tag           = "MAPPED_SKY_SCI_MXU";
00937             unmapped_sky_tag         = "UNMAPPED_SKY_SCI_MXU";
00938             object_table_tag         = "OBJECT_TABLE_SCI_MXU";
00939             reduced_sky_tag          = "REDUCED_SKY_SCI_MXU";
00940             reduced_error_tag        = "REDUCED_ERROR_SCI_MXU";
00941             reduced_flux_error_tag   = "REDUCED_FLUX_ERROR_SCI_MXU";
00942             specphot_tag             = "SPECPHOT_TABLE";
00943         }
00944 
00945         master_norm_flat_tag    = "MASTER_NORM_FLAT_MXU";
00946         disp_coeff_tag          = "DISP_COEFF_MXU";
00947         curv_coeff_tag          = "CURV_COEFF_MXU";
00948         slit_location_tag       = "SLIT_LOCATION_MXU";
00949         global_sky_spectrum_tag = "GLOBAL_SKY_SPECTRUM_MXU";
00950 
00951         if (!cpl_frameset_count_tags(frameset, master_norm_flat_tag)) {
00952             master_norm_flat_tag   = "MASTER_NORM_FLAT_LONG_MXU";
00953             disp_coeff_tag         = "DISP_COEFF_LONG_MXU";
00954             slit_location_tag      = "SLIT_LOCATION_LONG_MXU";
00955         }
00956     }
00957 
00958     if (lss) {
00959         cpl_msg_info(recipe, "LSS data found");
00960 
00961         if (cosmics && !skyglobal)
00962             fors_science_exit("Cosmic rays correction for LSS "
00963                               "data requires --skyglobal=true");
00964 
00965         if (standard) {
00966             science_tag              = "STANDARD_LSS";
00967             reduced_science_tag      = "REDUCED_STD_LSS";
00968             reduced_flux_science_tag = "REDUCED_FLUX_STD_LSS";
00969             unmapped_science_tag     = "UNMAPPED_STD_LSS";
00970             mapped_science_tag       = "MAPPED_STD_LSS";
00971             mapped_flux_science_tag  = "MAPPED_FLUX_STD_LSS";
00972             mapped_science_sky_tag   = "MAPPED_ALL_STD_LSS";
00973             skylines_offsets_tag     = "SKY_SHIFTS_LONG_STD_LSS";
00974             wavelength_map_sky_tag   = "WAVELENGTH_MAP_STD_LSS";
00975             disp_coeff_sky_tag       = "DISP_COEFF_STD_LSS";
00976             mapped_sky_tag           = "MAPPED_SKY_STD_LSS";
00977             unmapped_sky_tag         = "UNMAPPED_SKY_STD_LSS";
00978             object_table_tag         = "OBJECT_TABLE_STD_LSS";
00979             reduced_sky_tag          = "REDUCED_SKY_STD_LSS";
00980             reduced_error_tag        = "REDUCED_ERROR_STD_LSS";
00981             reduced_flux_error_tag   = "REDUCED_FLUX_ERROR_STD_LSS";
00982             specphot_tag             = "SPECPHOT_TABLE";
00983         }
00984         else {
00985             science_tag              = "SCIENCE_LSS";
00986             reduced_science_tag      = "REDUCED_SCI_LSS";
00987             reduced_flux_science_tag = "REDUCED_FLUX_SCI_LSS";
00988             unmapped_science_tag     = "UNMAPPED_SCI_LSS";
00989             mapped_science_tag       = "MAPPED_SCI_LSS";
00990             mapped_flux_science_tag  = "MAPPED_FLUX_SCI_LSS";
00991             mapped_science_sky_tag   = "MAPPED_ALL_SCI_LSS";
00992             skylines_offsets_tag     = "SKY_SHIFTS_LONG_SCI_LSS";
00993             wavelength_map_sky_tag   = "WAVELENGTH_MAP_SCI_LSS";
00994             disp_coeff_sky_tag       = "DISP_COEFF_SCI_LSS";
00995             mapped_sky_tag           = "MAPPED_SKY_SCI_LSS";
00996             unmapped_sky_tag         = "UNMAPPED_SKY_SCI_LSS";
00997             object_table_tag         = "OBJECT_TABLE_SCI_LSS";
00998             reduced_sky_tag          = "REDUCED_SKY_SCI_LSS";
00999             reduced_error_tag        = "REDUCED_ERROR_SCI_LSS";
01000             reduced_flux_error_tag   = "REDUCED_FLUX_ERROR_SCI_LSS";
01001             specphot_tag             = "SPECPHOT_TABLE";
01002         }
01003 
01004         master_norm_flat_tag    = "MASTER_NORM_FLAT_LSS";
01005         disp_coeff_tag          = "DISP_COEFF_LSS";
01006         slit_location_tag       = "SLIT_LOCATION_LSS";
01007         global_sky_spectrum_tag = "GLOBAL_SKY_SPECTRUM_LSS";
01008     }
01009 
01010     if (mos) {
01011         cpl_msg_info(recipe, "MOS data found");
01012         if (standard) {
01013             science_tag              = "STANDARD_MOS";
01014             reduced_science_tag      = "REDUCED_STD_MOS";
01015             reduced_flux_science_tag = "REDUCED_FLUX_STD_MOS";
01016             unmapped_science_tag     = "UNMAPPED_STD_MOS";
01017             mapped_science_tag       = "MAPPED_STD_MOS";
01018             mapped_flux_science_tag  = "MAPPED_FLUX_STD_MOS";
01019             mapped_science_sky_tag   = "MAPPED_ALL_STD_MOS";
01020             skylines_offsets_tag     = "SKY_SHIFTS_SLIT_STD_MOS";
01021             wavelength_map_sky_tag   = "WAVELENGTH_MAP_STD_MOS";
01022             disp_coeff_sky_tag       = "DISP_COEFF_STD_MOS";
01023             mapped_sky_tag           = "MAPPED_SKY_STD_MOS";
01024             unmapped_sky_tag         = "UNMAPPED_SKY_STD_MOS";
01025             object_table_tag         = "OBJECT_TABLE_STD_MOS";
01026             reduced_sky_tag          = "REDUCED_SKY_STD_MOS";
01027             reduced_error_tag        = "REDUCED_ERROR_STD_MOS";
01028             reduced_flux_error_tag   = "REDUCED_FLUX_ERROR_STD_MOS";
01029             specphot_tag             = "SPECPHOT_TABLE";
01030         }
01031         else {
01032             science_tag              = "SCIENCE_MOS";
01033             reduced_science_tag      = "REDUCED_SCI_MOS";
01034             reduced_flux_science_tag = "REDUCED_FLUX_SCI_MOS";
01035             unmapped_science_tag     = "UNMAPPED_SCI_MOS";
01036             mapped_science_tag       = "MAPPED_SCI_MOS";
01037             mapped_flux_science_tag  = "MAPPED_FLUX_SCI_MOS";
01038             mapped_science_sky_tag   = "MAPPED_ALL_SCI_MOS";
01039             skylines_offsets_tag     = "SKY_SHIFTS_SLIT_SCI_MOS";
01040             wavelength_map_sky_tag   = "WAVELENGTH_MAP_SCI_MOS";
01041             disp_coeff_sky_tag       = "DISP_COEFF_SCI_MOS";
01042             mapped_sky_tag           = "MAPPED_SKY_SCI_MOS";
01043             unmapped_sky_tag         = "UNMAPPED_SKY_SCI_MOS";
01044             object_table_tag         = "OBJECT_TABLE_SCI_MOS";
01045             reduced_sky_tag          = "REDUCED_SKY_SCI_MOS";
01046             reduced_error_tag        = "REDUCED_ERROR_SCI_MOS";
01047             reduced_flux_error_tag   = "REDUCED_FLUX_ERROR_SCI_MOS";
01048             specphot_tag             = "SPECPHOT_TABLE";
01049         }
01050 
01051         master_norm_flat_tag    = "MASTER_NORM_FLAT_MOS";
01052         disp_coeff_tag          = "DISP_COEFF_MOS";
01053         curv_coeff_tag          = "CURV_COEFF_MOS";
01054         slit_location_tag       = "SLIT_LOCATION_MOS";
01055         global_sky_spectrum_tag = "GLOBAL_SKY_SPECTRUM_MOS";
01056 
01057         if (!cpl_frameset_count_tags(frameset, master_norm_flat_tag)) {
01058             master_norm_flat_tag   = "MASTER_NORM_FLAT_LONG_MOS";
01059             disp_coeff_tag         = "DISP_COEFF_LONG_MOS";
01060             slit_location_tag      = "SLIT_LOCATION_LONG_MOS";
01061         }
01062     }
01063 
01064     if (cpl_frameset_count_tags(frameset, "MASTER_BIAS") == 0)
01065         fors_science_exit("Missing required input: MASTER_BIAS");
01066 
01067     if (cpl_frameset_count_tags(frameset, "MASTER_BIAS") > 1)
01068         fors_science_exit("Too many in input: MASTER_BIAS");
01069 
01070     if (skyalign >= 0)
01071         if (cpl_frameset_count_tags(frameset, "MASTER_SKYLINECAT") > 1)
01072             fors_science_exit("Too many in input: MASTER_SKYLINECAT");
01073 
01074     if (cpl_frameset_count_tags(frameset, disp_coeff_tag) == 0) {
01075         cpl_msg_error(recipe, "Missing required input: %s", disp_coeff_tag);
01076         fors_science_exit(NULL);
01077     }
01078 
01079     if (cpl_frameset_count_tags(frameset, disp_coeff_tag) > 1) {
01080         cpl_msg_error(recipe, "Too many in input: %s", disp_coeff_tag);
01081         fors_science_exit(NULL);
01082     }
01083 
01084     if (cpl_frameset_count_tags(frameset, slit_location_tag) == 0) {
01085         cpl_msg_error(recipe, "Missing required input: %s",
01086                       slit_location_tag);
01087         fors_science_exit(NULL);
01088     }
01089 
01090     if (cpl_frameset_count_tags(frameset, slit_location_tag) > 1) {
01091         cpl_msg_error(recipe, "Too many in input: %s", slit_location_tag);
01092         fors_science_exit(NULL);
01093     }
01094 
01095     if (cpl_frameset_count_tags(frameset, master_norm_flat_tag) > 1) {
01096         if (flatfield) {
01097             cpl_msg_error(recipe, "Too many in input: %s", 
01098                           master_norm_flat_tag);
01099             fors_science_exit(NULL);
01100         }
01101         else {
01102             cpl_msg_warning(recipe, "%s in input are ignored, "
01103                             "since flat field correction was not requested", 
01104                             master_norm_flat_tag);
01105         }
01106     }
01107 
01108     if (cpl_frameset_count_tags(frameset, master_norm_flat_tag) == 1) {
01109         if (!flatfield) {
01110             cpl_msg_warning(recipe, "%s in input is ignored, "
01111                             "since flat field correction was not requested", 
01112                             master_norm_flat_tag);
01113         }
01114     }
01115 
01116     if (cpl_frameset_count_tags(frameset, master_norm_flat_tag) == 0) {
01117         if (flatfield) {
01118             cpl_msg_error(recipe, "Flat field correction was requested, "
01119                           "but no %s are found in input",
01120                           master_norm_flat_tag);
01121             fors_science_exit(NULL);
01122         }
01123     }
01124 
01125     if (standard) {
01126     
01127         if (cpl_frameset_count_tags(frameset, "EXTINCT_TABLE") == 0) {
01128             cpl_msg_warning(recipe, "An EXTINCT_TABLE was not found in input: "
01129                             "instrument response curve will not be produced.");
01130             standard = 0;
01131         }
01132 
01133         if (cpl_frameset_count_tags(frameset, "EXTINCT_TABLE") > 1)
01134             fors_science_exit("Too many in input: EXTINCT_TABLE");
01135     
01136         if (cpl_frameset_count_tags(frameset, "STD_FLUX_TABLE") == 0) {
01137             cpl_msg_warning(recipe, "A STD_FLUX_TABLE was not found in input: "
01138                             "instrument response curve will not be produced.");
01139             standard = 0;
01140         }
01141 
01142         if (cpl_frameset_count_tags(frameset, "STD_FLUX_TABLE") > 1)
01143             fors_science_exit("Too many in input: STD_FLUX_TABLE");
01144 
01145         if (!dfs_equal_keyword(frameset, "ESO OBS TARG NAME")) {
01146             cpl_msg_warning(recipe, "The target name of observation does not "
01147                             "match the standard star catalog: "
01148                             "instrument response curve will not be produced.");
01149             standard = 0;
01150         }
01151     }
01152 
01153     if (photometry) {
01154         if (cpl_frameset_count_tags(frameset, "EXTINCT_TABLE") == 0) {
01155             fors_science_exit("An EXTINCT_TABLE was not found in input: "
01156                               "the requested photometric calibrated spectra "
01157                               "cannot be produced.");
01158         }
01159 
01160         if (cpl_frameset_count_tags(frameset, "EXTINCT_TABLE") > 1)
01161             fors_science_exit("Too many in input: EXTINCT_TABLE");
01162 
01163         have_phot  = cpl_frameset_count_tags(frameset, specphot_tag);
01164         have_phot += cpl_frameset_count_tags(frameset, master_specphot_tag);
01165 
01166         if (!standard) {
01167             if (have_phot == 0) {
01168                 cpl_msg_warning(recipe, 
01169                                 "A SPECPHOT_TABLE was not found in input: "
01170                                 "the requested photometric calibrated "
01171                                 "spectra cannot be produced.");
01172                 photometry = 0;
01173             }
01174 
01175             if (have_phot > 1)
01176                 fors_science_exit("Too many in input: SPECPHOT_TABLE");
01177         }
01178     }
01179 
01180     cpl_msg_indent_less();
01181 
01182 
01183     /*
01184      * Loading input data
01185      */
01186 
01187     exptime = cpl_calloc(nscience, sizeof(double));
01188 
01189     if (nscience > 1) {
01190 
01191         cpl_msg_info(recipe, "Load %d scientific frames and median them...",
01192                      nscience);
01193         cpl_msg_indent_more();
01194 
01195         all_science = cpl_imagelist_new();
01196 
01197         header = dfs_load_header(frameset, science_tag, 0);
01198 
01199         if (header == NULL)
01200             fors_science_exit("Cannot load scientific frame header");
01201 
01202         alltime = exptime[0] = cpl_propertylist_get_double(header, "EXPTIME");
01203 
01204         if (cpl_error_get_code() != CPL_ERROR_NONE)
01205             fors_science_exit("Missing keyword EXPTIME in scientific "
01206                               "frame header");
01207 
01208         if (standard || photometry) {
01209             airmass = fors_get_airmass(header);
01210             if (airmass < 0.0) 
01211                 fors_science_exit("Missing airmass information in "
01212                                   "scientific frame header");
01213         }
01214 
01215         cpl_propertylist_delete(header); header = NULL;
01216 
01217         cpl_msg_info(recipe, "Scientific frame 1 exposure time: %.2f s", 
01218                      exptime[0]);
01219 
01220         for (i = 1; i < nscience; i++) {
01221 
01222             header = dfs_load_header(frameset, NULL, 0);
01223 
01224             if (header == NULL)
01225                 fors_science_exit("Cannot load scientific frame header");
01226     
01227             exptime[i] = cpl_propertylist_get_double(header, "EXPTIME");
01228 
01229             alltime += exptime[i];
01230     
01231             if (cpl_error_get_code() != CPL_ERROR_NONE)
01232                 fors_science_exit("Missing keyword EXPTIME in scientific "
01233                                   "frame header");
01234     
01235             cpl_propertylist_delete(header); header = NULL;
01236 
01237             cpl_msg_info(recipe, "Scientific frame %d exposure time: %.2f s", 
01238                          i+1, exptime[i]);
01239         }
01240 
01241         spectra = dfs_load_image(frameset, science_tag, CPL_TYPE_FLOAT, 0, 0);
01242 
01243         if (spectra == NULL)
01244             fors_science_exit("Cannot load scientific frame");
01245 
01246         cpl_image_divide_scalar(spectra, exptime[0]);
01247         cpl_imagelist_set(all_science, spectra, 0); spectra = NULL;
01248 
01249         for (i = 1; i < nscience; i++) {
01250 
01251             spectra = dfs_load_image(frameset, NULL, CPL_TYPE_FLOAT, 0, 0);
01252 
01253             if (spectra) {
01254                 cpl_image_divide_scalar(spectra, exptime[i]);
01255                 cpl_imagelist_set(all_science, spectra, i); spectra = NULL;
01256             }
01257             else
01258                 fors_science_exit("Cannot load scientific frame");
01259 
01260         }
01261 
01262         spectra = cpl_imagelist_collapse_median_create(all_science);
01263         cpl_image_multiply_scalar(spectra, alltime);
01264 
01265         cpl_imagelist_delete(all_science);
01266     }
01267     else {
01268         cpl_msg_info(recipe, "Load scientific exposure...");
01269         cpl_msg_indent_more();
01270 
01271         header = dfs_load_header(frameset, science_tag, 0);
01272 
01273         if (header == NULL)
01274             fors_science_exit("Cannot load scientific frame header");
01275 
01276         if (standard || photometry) {
01277             airmass = fors_get_airmass(header);
01278             if (airmass < 0.0) 
01279                 fors_science_exit("Missing airmass information in "
01280                                   "scientific frame header");
01281         }
01282 
01283         /*
01284          * Insert here a check on supported filters:
01285          */
01286 
01287         wheel4 = (char *)cpl_propertylist_get_string(header, 
01288                                                      "ESO INS OPTI9 TYPE");
01289         if (cpl_error_get_code() != CPL_ERROR_NONE) {
01290             fors_science_exit("Missing ESO INS OPTI9 TYPE in flat header");
01291         }
01292 
01293         if (strcmp("FILT", wheel4) == 0) {
01294             wheel4 = (char *)cpl_propertylist_get_string(header,
01295                                                          "ESO INS OPTI9 NAME");
01296             cpl_msg_error(recipe, "Unsupported filter: %s", wheel4);
01297             fors_science_exit(NULL);
01298         }
01299 
01300         alltime = exptime[0] = cpl_propertylist_get_double(header, "EXPTIME");
01301 
01302         if (cpl_error_get_code() != CPL_ERROR_NONE)
01303             fors_science_exit("Missing keyword EXPTIME in scientific "
01304                               "frame header");
01305 
01306         cpl_propertylist_delete(header); header = NULL;
01307 
01308         cpl_msg_info(recipe, "Scientific frame exposure time: %.2f s", 
01309                      exptime[0]);
01310 
01311         spectra = dfs_load_image(frameset, science_tag, CPL_TYPE_FLOAT, 0, 0);
01312     }
01313 
01314     if (spectra == NULL)
01315         fors_science_exit("Cannot load scientific frame");
01316 
01317     cpl_free(exptime); exptime = NULL;
01318 
01319     cpl_msg_indent_less();
01320 
01321 
01322     /*
01323      * Get the reference wavelength and the rebin factor along the
01324      * dispersion direction from a scientific exposure
01325      */
01326 
01327     header = dfs_load_header(frameset, science_tag, 0);
01328 
01329     if (header == NULL)
01330         fors_science_exit("Cannot load scientific frame header");
01331 
01332     instrume = (char *)cpl_propertylist_get_string(header, "INSTRUME");
01333     if (instrume == NULL)
01334         fors_science_exit("Missing keyword INSTRUME in scientific header");
01335     instrume = cpl_strdup(instrume);
01336 
01337     if (instrume[4] == '1')
01338         snprintf(version, 80, "%s/%s", "fors1", VERSION);
01339     if (instrume[4] == '2')
01340         snprintf(version, 80, "%s/%s", "fors2", VERSION);
01341 
01342     reference = cpl_propertylist_get_double(header, "ESO INS GRIS1 WLEN");
01343 
01344     if (cpl_error_get_code() != CPL_ERROR_NONE)
01345         fors_science_exit("Missing keyword ESO INS GRIS1 WLEN in scientific "
01346                         "frame header");
01347 
01348     if (reference < 3000.0)   /* Perhaps in nanometers... */
01349         reference *= 10;
01350 
01351     if (reference < 3000.0 || reference > 13000.0) {
01352         cpl_msg_error(recipe, "Invalid central wavelength %.2f read from "
01353                       "keyword ESO INS GRIS1 WLEN in scientific frame header",
01354                       reference);
01355         fors_science_exit(NULL);
01356     }
01357 
01358     cpl_msg_info(recipe, "The central wavelength is: %.2f", reference);
01359 
01360     rebin = cpl_propertylist_get_int(header, "ESO DET WIN1 BINX");
01361 
01362     if (cpl_error_get_code() != CPL_ERROR_NONE)
01363         fors_science_exit("Missing keyword ESO DET WIN1 BINX in scientific "
01364                         "frame header");
01365 
01366     if (rebin != 1) {
01367         dispersion *= rebin;
01368         cpl_msg_warning(recipe, "The rebin factor is %d, and therefore the "
01369                         "resampling step used is %f A/pixel", rebin, 
01370                         dispersion);
01371         ext_radius /= rebin;
01372         cpl_msg_warning(recipe, "The rebin factor is %d, and therefore the "
01373                         "extraction radius used is %d pixel", rebin, 
01374                         ext_radius);
01375     }
01376 
01377     gain = cpl_propertylist_get_double(header, "ESO DET OUT1 CONAD");
01378 
01379     if (cpl_error_get_code() != CPL_ERROR_NONE)
01380         fors_science_exit("Missing keyword ESO DET OUT1 CONAD in scientific "
01381                           "frame header");
01382 
01383     cpl_msg_info(recipe, "The gain factor is: %.2f e-/ADU", gain);
01384 
01385     ron = cpl_propertylist_get_double(header, "ESO DET OUT1 RON");
01386 
01387     if (cpl_error_get_code() != CPL_ERROR_NONE)
01388         fors_science_exit("Missing keyword ESO DET OUT1 RON in scientific "
01389                           "frame header");
01390 
01391     ron /= gain;     /* Convert from electrons to ADU */
01392 
01393     cpl_msg_info(recipe, "The read-out-noise is: %.2f ADU", ron);
01394 
01395     if (mos || mxu) {
01396 /* goto skip; */
01397         if (mos)
01398             maskslits = mos_load_slits_fors_mos(header);
01399         else
01400             maskslits = mos_load_slits_fors_mxu(header);
01401 
01402         /*
01403          * Check if all slits have the same X offset: in such case,
01404          * treat the observation as a long-slit one!
01405          */
01406 
01407         mxpos = cpl_table_get_column_median(maskslits, "xtop");
01408         xpos = cpl_table_get_data_double(maskslits, "xtop");
01409         nslits = cpl_table_get_nrow(maskslits);
01410      
01411         treat_as_lss = 1;
01412         for (i = 0; i < nslits; i++) { 
01413             if (fabs(mxpos-xpos[i]) > 0.01) {
01414                 treat_as_lss = 0;
01415                 break;
01416             }
01417         }
01418 
01419         cpl_table_delete(maskslits); maskslits = NULL;
01420 /* skip: */
01421 
01422         if (treat_as_lss) {
01423             cpl_msg_warning(recipe, "All MOS slits have the same offset: %.2f\n"
01424                             "The LSS data reduction strategy is applied!",
01425                             mxpos);
01426             if (mos) {
01427                 if (standard) {
01428                     skylines_offsets_tag   = "SKY_SHIFTS_LONG_STD_MOS";
01429                 }
01430                 else {
01431                     skylines_offsets_tag   = "SKY_SHIFTS_LONG_SCI_MOS";
01432                 }
01433             }
01434             else {
01435                 if (standard) {
01436                     skylines_offsets_tag   = "SKY_SHIFTS_LONG_STD_MXU";
01437                 }
01438                 else {
01439                     skylines_offsets_tag   = "SKY_SHIFTS_LONG_SCI_MXU";
01440                 }
01441             }
01442         }
01443     }
01444 
01445     if (lss || treat_as_lss) {
01446         if (skylocal) {
01447             if (cosmics)
01448                 fors_science_exit("Cosmic rays correction for LSS or LSS-like "
01449                                   "data requires --skyglobal=true");
01450             skymedian = skylocal;
01451             skylocal = 0;
01452         }
01453     }
01454     else {
01455         if (cpl_frameset_count_tags(frameset, curv_coeff_tag) == 0) {
01456             cpl_msg_error(recipe, "Missing required input: %s", curv_coeff_tag);
01457             fors_science_exit(NULL);
01458         }
01459 
01460         if (cpl_frameset_count_tags(frameset, curv_coeff_tag) > 1) {
01461             cpl_msg_error(recipe, "Too many in input: %s", curv_coeff_tag);
01462             fors_science_exit(NULL);
01463         }
01464     }
01465 
01466     /* Leave the header on for the next step... */
01467 
01468 
01469     /*
01470      * Remove the master bias
01471      */
01472 
01473     cpl_msg_info(recipe, "Remove the master bias...");
01474 
01475     bias = dfs_load_image(frameset, "MASTER_BIAS", CPL_TYPE_FLOAT, 0, 1);
01476 
01477     if (bias == NULL)
01478         fors_science_exit("Cannot load master bias");
01479 
01480     overscans = mos_load_overscans_vimos(header, 1);
01481     cpl_propertylist_delete(header); header = NULL;
01482     dummy = mos_remove_bias(spectra, bias, overscans);
01483     cpl_image_delete(spectra); spectra = dummy; dummy = NULL;
01484     cpl_image_delete(bias); bias = NULL;
01485     cpl_table_delete(overscans); overscans = NULL;
01486 
01487     if (spectra == NULL)
01488         fors_science_exit("Cannot remove bias from scientific frame");
01489 
01490     ccd_xsize = nx = cpl_image_get_size_x(spectra);
01491     ccd_ysize = ny = cpl_image_get_size_y(spectra);
01492 
01493     cpl_msg_indent_less();
01494     cpl_msg_info(recipe, "Load normalised flat field (if present)...");
01495     cpl_msg_indent_more();
01496 
01497     if (flatfield) {
01498 
01499         norm_flat = dfs_load_image(frameset, master_norm_flat_tag, 
01500                                    CPL_TYPE_FLOAT, 0, 1);
01501 
01502         if (norm_flat) {
01503             cpl_msg_info(recipe, "Apply flat field correction...");
01504             if (cpl_image_divide(spectra, norm_flat) != CPL_ERROR_NONE) {
01505                 cpl_msg_error(recipe, "Failure of flat field correction: %s",
01506                               cpl_error_get_message());
01507                 fors_science_exit(NULL);
01508             }
01509             cpl_image_delete(norm_flat); norm_flat = NULL;
01510         }
01511         else {
01512             cpl_msg_error(recipe, "Cannot load input %s for flat field "
01513                           "correction", master_norm_flat_tag);
01514             fors_science_exit(NULL);
01515         }
01516 
01517     }
01518 
01519 
01520     if (skyalign >= 0) {
01521         cpl_msg_indent_less();
01522         cpl_msg_info(recipe, "Load input sky line catalog...");
01523         cpl_msg_indent_more();
01524 
01525         wavelengths = dfs_load_table(frameset, "MASTER_SKYLINECAT", 1);
01526 
01527         if (wavelengths) {
01528 
01529             /*
01530              * Cast the wavelengths into a (double precision) CPL vector
01531              */
01532 
01533             nlines = cpl_table_get_nrow(wavelengths);
01534 
01535             if (nlines == 0)
01536                 fors_science_exit("Empty input sky line catalog");
01537 
01538             if (cpl_table_has_column(wavelengths, wcolumn) != 1) {
01539                 cpl_msg_error(recipe, "Missing column %s in input line "
01540                               "catalog table", wcolumn);
01541                 fors_science_exit(NULL);
01542             }
01543 
01544             line = cpl_malloc(nlines * sizeof(double));
01545     
01546             for (i = 0; i < nlines; i++)
01547                 line[i] = cpl_table_get(wavelengths, wcolumn, i, NULL);
01548 
01549             cpl_table_delete(wavelengths); wavelengths = NULL;
01550 
01551             lines = cpl_vector_wrap(nlines, line);
01552         }
01553         else {
01554             cpl_msg_info(recipe, "No sky line catalog found in input - fine!");
01555         }
01556     }
01557 
01558 
01559     /*
01560      * Load the spectral curvature table, or provide a dummy one
01561      * in case of LSS or LSS-like data (single slit with flat curvature)
01562      */
01563 
01564     if (!(lss || treat_as_lss)) {
01565         polytraces = dfs_load_table(frameset, curv_coeff_tag, 1);
01566         if (polytraces == NULL)
01567             fors_science_exit("Cannot load spectral curvature table");
01568     }
01569 
01570 
01571     /*
01572      * Load the slit location table
01573      */
01574 
01575     slits = dfs_load_table(frameset, slit_location_tag, 1);
01576     if (slits == NULL)
01577         fors_science_exit("Cannot load slits location table");
01578 
01579     if (lss || treat_as_lss) {
01580         int first_row = cpl_table_get_double(slits, "ybottom", 0, NULL);
01581         int last_row = cpl_table_get_double(slits, "ytop", 0, NULL);
01582         int ylow, yhig;
01583 
01584         ylow = first_row + 1;
01585         yhig = last_row + 1;
01586 
01587         dummy = cpl_image_extract(spectra, 1, ylow, nx, yhig);
01588         cpl_image_delete(spectra); spectra = dummy; dummy = NULL;
01589         ny = cpl_image_get_size_y(spectra);
01590     }
01591 
01592 
01593     /*
01594      * Load the wavelength calibration table
01595      */
01596 
01597     idscoeff = dfs_load_table(frameset, disp_coeff_tag, 1);
01598     if (idscoeff == NULL)
01599         fors_science_exit("Cannot load wavelength calibration table");
01600 
01601     cpl_msg_indent_less();
01602     cpl_msg_info(recipe, "Processing scientific spectra...");
01603     cpl_msg_indent_more();
01604 
01605     /*
01606      * This one will also generate the spatial map from the spectral 
01607      * curvature table (in the case of multislit data)
01608      */
01609 
01610     if (lss || treat_as_lss) {
01611         smapped = cpl_image_duplicate(spectra);
01612     }
01613     else {
01614         coordinate = cpl_image_new(nx, ny, CPL_TYPE_FLOAT);
01615 
01616         smapped = mos_spatial_calibration(spectra, slits, polytraces, reference,
01617                                           startwavelength, endwavelength,
01618                                           dispersion, flux, coordinate);
01619     }
01620 
01621 
01622     /*
01623      * Generate a rectified wavelength map from the wavelength calibration 
01624      * table
01625      */
01626 
01627     rainbow = mos_map_idscoeff(idscoeff, nx, reference, startwavelength, 
01628                                endwavelength);
01629 
01630     if (dispersion > 1.0)
01631         highres = 0;
01632     else
01633         highres = 1;
01634 
01635     if (skyalign >= 0) {
01636         if (skyalign) {
01637             cpl_msg_info(recipe, "Align wavelength solution to reference "
01638             "skylines applying %d order residual fit...", skyalign);
01639         }
01640         else {
01641             cpl_msg_info(recipe, "Align wavelength solution to reference "
01642             "skylines applying median offset...");
01643         }
01644 
01645         if (lss || treat_as_lss) {
01646             offsets = mos_wavelength_align_lss(smapped, reference, 
01647                                                startwavelength, endwavelength, 
01648                                                idscoeff, lines, highres, 
01649                                                skyalign, rainbow, 4);
01650         }
01651         else {
01652             offsets = mos_wavelength_align(smapped, slits, reference, 
01653                                            startwavelength, endwavelength, 
01654                                            idscoeff, lines, highres, skyalign, 
01655                                            rainbow, 4);
01656         }
01657 
01658         if (offsets) {
01659             if (standard)
01660                 cpl_msg_warning(recipe, "Alignment of the wavelength solution "
01661                                 "to reference sky lines may be unreliable in "
01662                                 "this case!");
01663 
01664             if (dfs_save_table(frameset, offsets, skylines_offsets_tag, NULL, 
01665                                parlist, recipe, version))
01666                 fors_science_exit(NULL);
01667 
01668             cpl_table_delete(offsets); offsets = NULL;
01669         }
01670         else {
01671             if (cpl_error_get_code()) {
01672                 if (cpl_error_get_code() == CPL_ERROR_INCOMPATIBLE_INPUT) {
01673                     cpl_msg_error(recipe, "The IDS coeff table is "
01674                     "incompatible with the input slit position table.");
01675                 }
01676                 cpl_msg_error(cpl_error_get_where(), cpl_error_get_message());
01677                 fors_science_exit(NULL);
01678             }
01679             cpl_msg_warning(recipe, "Alignment of the wavelength solution "
01680                             "to reference sky lines could not be done!");
01681             skyalign = -1;
01682         }
01683 
01684     }
01685 
01686     if (lss || treat_as_lss) {
01687         int first_row = cpl_table_get_double(slits, "ybottom", 0, NULL);
01688         int last_row = cpl_table_get_double(slits, "ytop", 0, NULL);
01689         int ylow, yhig;
01690 
01691         ylow = first_row + 1;
01692         yhig = last_row + 1;
01693 
01694         wavemap = cpl_image_new(ccd_xsize, ccd_ysize, CPL_TYPE_FLOAT);
01695         cpl_image_copy(wavemap, rainbow, 1, ylow);
01696 
01697         wavemaplss = cpl_image_extract(wavemap, 1, ylow, nx, yhig);
01698     }
01699     else {
01700         wavemap = mos_map_wavelengths(coordinate, rainbow, slits, 
01701                                       polytraces, reference, 
01702                                       startwavelength, endwavelength,
01703                                       dispersion);
01704     }
01705 
01706     cpl_image_delete(rainbow); rainbow = NULL;
01707     cpl_image_delete(coordinate); coordinate = NULL;
01708 
01709     /*
01710      * Here the wavelength calibrated slit spectra are created. This frame
01711      * contains sky_science.
01712      */
01713 
01714     mapped_sky = mos_wavelength_calibration(smapped, reference,
01715                                             startwavelength, endwavelength,
01716                                             dispersion, idscoeff, flux);
01717 
01718     cpl_msg_indent_less();
01719     cpl_msg_info(recipe, "Check applied wavelength against skylines...");
01720     cpl_msg_indent_more();
01721 
01722     mean_rms = mos_distortions_rms(mapped_sky, lines, startwavelength,
01723                                    dispersion, 6, highres);
01724 
01725     cpl_vector_delete(lines); lines = NULL;
01726 
01727     cpl_msg_info(recipe, "Mean residual: %f", mean_rms);
01728 
01729     mean_rms = cpl_table_get_column_mean(idscoeff, "error");
01730 
01731     cpl_msg_info(recipe, "Mean model accuracy: %f pixel (%f A)",
01732                  mean_rms, mean_rms * dispersion);
01733 
01734     header = cpl_propertylist_new();
01735     cpl_propertylist_update_double(header, "CRPIX1", 1.0);
01736     cpl_propertylist_update_double(header, "CRPIX2", 1.0);
01737     cpl_propertylist_update_double(header, "CRVAL1", 
01738                                    startwavelength + dispersion/2);
01739     cpl_propertylist_update_double(header, "CRVAL2", 1.0);
01740     /* cpl_propertylist_update_double(header, "CDELT1", dispersion);
01741     cpl_propertylist_update_double(header, "CDELT2", 1.0); */
01742     cpl_propertylist_update_double(header, "CD1_1", dispersion);
01743     cpl_propertylist_update_double(header, "CD1_2", 0.0);
01744     cpl_propertylist_update_double(header, "CD2_1", 0.0);
01745     cpl_propertylist_update_double(header, "CD2_2", 1.0);
01746     cpl_propertylist_update_string(header, "CTYPE1", "LINEAR");
01747     cpl_propertylist_update_string(header, "CTYPE2", "PIXEL");
01748 
01749     if (time_normalise) {
01750         cpl_propertylist_update_string(header, "BUNIT", "ADU/s");
01751         dummy = cpl_image_divide_scalar_create(mapped_sky, alltime);
01752         if (dfs_save_image(frameset, dummy, mapped_science_sky_tag, header, 
01753                            parlist, recipe, version))
01754             fors_science_exit(NULL);
01755         cpl_image_delete(dummy); dummy = NULL;
01756     }
01757     else {
01758         cpl_propertylist_update_string(header, "BUNIT", "ADU");
01759         if (dfs_save_image(frameset, mapped_sky, mapped_science_sky_tag, 
01760                            header, parlist, recipe, version))
01761             fors_science_exit(NULL);
01762     }
01763 
01764 /*    if (skyglobal == 0 && skymedian < 0) {    NSS */
01765     if (skyglobal == 0 && skymedian == 0 && skylocal == 0) {
01766         cpl_image_delete(mapped_sky); mapped_sky = NULL;
01767     }
01768 
01769     if (skyglobal || skylocal) {
01770 
01771         cpl_msg_indent_less();
01772 
01773         if (skyglobal) {
01774             cpl_msg_info(recipe, "Global sky determination...");
01775             cpl_msg_indent_more();
01776             skymap = cpl_image_new(nx, ny, CPL_TYPE_FLOAT);
01777             if (lss || treat_as_lss) {
01778                 sky = mos_sky_map_super(spectra, wavemaplss, dispersion, 
01779                                         2.0, 50, skymap);
01780             }
01781             else {
01782                 sky = mos_sky_map_super(spectra, wavemap, dispersion, 
01783                                         2.0, 50, skymap);
01784             }
01785             if (sky) {
01786                 cpl_image_subtract(spectra, skymap);
01787             }
01788             else {
01789                 cpl_image_delete(skymap); skymap = NULL;
01790             }
01791         }
01792         else {
01793             cpl_msg_info(recipe, "Local sky determination...");
01794             cpl_msg_indent_more();
01795             skymap = mos_subtract_sky(spectra, slits, polytraces, reference,
01796                            startwavelength, endwavelength, dispersion);
01797         }
01798 
01799         if (skymap) {
01800             if (skyglobal) {
01801                 if (time_normalise)
01802                     cpl_table_divide_scalar(sky, "sky", alltime);
01803                 if (dfs_save_table(frameset, sky, global_sky_spectrum_tag, 
01804                                    NULL, parlist, recipe, version))
01805                     fors_science_exit(NULL);
01806     
01807                 cpl_table_delete(sky); sky = NULL;
01808             }
01809 
01810             save_header = dfs_load_header(frameset, science_tag, 0);
01811 
01812             if (time_normalise)
01813                 cpl_image_divide_scalar(skymap, alltime);
01814             if (dfs_save_image(frameset, skymap, unmapped_sky_tag,
01815                                save_header, parlist, recipe, version))
01816                 fors_science_exit(NULL);
01817 
01818             cpl_image_delete(skymap); skymap = NULL;
01819 
01820             if (dfs_save_image(frameset, spectra, unmapped_science_tag,
01821                                save_header, parlist, recipe, version))
01822                 fors_science_exit(NULL);
01823 
01824             cpl_propertylist_delete(save_header); save_header = NULL;
01825 
01826             if (cosmics) {
01827                 cpl_msg_info(recipe, "Removing cosmic rays...");
01828                 mos_clean_cosmics(spectra, gain, -1., -1.);
01829             }
01830 
01831             /*
01832              * The spatially rectified image, that contained the sky,
01833              * is replaced by a sky-subtracted spatially rectified image:
01834              */
01835 
01836             cpl_image_delete(smapped); smapped = NULL;
01837 
01838             if (lss || treat_as_lss) {
01839                 smapped = cpl_image_duplicate(spectra);
01840             }
01841             else {
01842                 smapped = mos_spatial_calibration(spectra, slits, polytraces, 
01843                                                   reference, startwavelength, 
01844                                                   endwavelength, dispersion, 
01845                                                   flux, NULL);
01846             }
01847         }
01848         else {
01849             cpl_msg_warning(recipe, "Sky subtraction failure");
01850             if (cosmics)
01851                 cpl_msg_warning(recipe, "Cosmic rays removal not performed!");
01852             cosmics = skylocal = skyglobal = 0;
01853         }
01854     }
01855 
01856     cpl_image_delete(spectra); spectra = NULL;
01857     cpl_table_delete(polytraces); polytraces = NULL;
01858     if (lss || treat_as_lss) 
01859         cpl_image_delete(wavemaplss); wavemaplss = NULL;
01860 
01861     if (skyalign >= 0) {
01862         save_header = dfs_load_header(frameset, science_tag, 0);
01863         cpl_propertylist_update_string(save_header, "BUNIT", "Angstrom");
01864         if (dfs_save_image(frameset, wavemap, wavelength_map_sky_tag,
01865                            save_header, parlist, recipe, version))
01866             fors_science_exit(NULL);
01867         cpl_propertylist_delete(save_header); save_header = NULL;
01868     }
01869 
01870     cpl_image_delete(wavemap); wavemap = NULL;
01871 
01872     mapped = mos_wavelength_calibration(smapped, reference,
01873                                         startwavelength, endwavelength,
01874                                         dispersion, idscoeff, flux);
01875 
01876     cpl_image_delete(smapped); smapped = NULL;
01877 
01878 /*    if (skymedian >= 0) {    NSS */
01879     if (skymedian) {
01880             cpl_msg_indent_less();
01881             cpl_msg_info(recipe, "Local sky determination...");
01882             cpl_msg_indent_more();
01883        
01884 /*   NSS      skylocalmap = mos_sky_local(mapped, slits, skymedian); */
01885 /*            skylocalmap = mos_sky_local(mapped, slits, 0);        */
01886             skylocalmap = mos_sky_local_old(mapped, slits);       
01887             cpl_image_subtract(mapped, skylocalmap);
01888 /*
01889             if (dfs_save_image(frameset, skylocalmap, mapped_sky_tag, header, 
01890                                parlist, recipe, version))
01891                 fors_science_exit(NULL);
01892 */
01893             cpl_image_delete(skylocalmap); skylocalmap = NULL;
01894     }
01895 
01896 /*    if (skyglobal || skymedian >= 0 || skylocal) {   NSS */
01897     if (skyglobal || skymedian || skylocal) {
01898 
01899         skylocalmap = cpl_image_subtract_create(mapped_sky, mapped);
01900 
01901         cpl_image_delete(mapped_sky); mapped_sky = NULL;
01902 
01903         if (time_normalise) {
01904             cpl_propertylist_update_string(header, "BUNIT", "ADU/s");
01905             dummy = cpl_image_divide_scalar_create(skylocalmap, alltime);
01906             if (dfs_save_image(frameset, dummy, mapped_sky_tag, header,
01907                                parlist, recipe, version))
01908                 fors_science_exit(NULL);
01909             cpl_image_delete(dummy); dummy = NULL;
01910         }
01911         else {
01912             cpl_propertylist_update_string(header, "BUNIT", "ADU");
01913             if (dfs_save_image(frameset, skylocalmap, mapped_sky_tag, header,
01914                                parlist, recipe, version))
01915                 fors_science_exit(NULL);
01916         }
01917 
01918         cpl_msg_indent_less();
01919         cpl_msg_info(recipe, "Object detection...");
01920         cpl_msg_indent_more();
01921 
01922         if (cosmics || nscience > 1) {
01923             dummy = mos_detect_objects(mapped, slits, slit_margin, ext_radius, 
01924                                        cont_radius);
01925         }
01926         else {
01927             mapped_cleaned = cpl_image_duplicate(mapped);
01928             mos_clean_cosmics(mapped_cleaned, gain, -1., -1.);
01929             dummy = mos_detect_objects(mapped_cleaned, slits, slit_margin, 
01930                                        ext_radius, cont_radius);
01931 
01932             cpl_image_delete(mapped_cleaned); mapped_cleaned = NULL;
01933         }
01934 
01935         cpl_image_delete(dummy); dummy = NULL;
01936 
01937         if (dfs_save_table(frameset, slits, object_table_tag, NULL, parlist, 
01938                            recipe, version))
01939             fors_science_exit(NULL);
01940 
01941         cpl_msg_indent_less();
01942         cpl_msg_info(recipe, "Object extraction...");
01943         cpl_msg_indent_more();
01944 
01945         images = mos_extract_objects(mapped, skylocalmap, slits, 
01946                                      ext_mode, ron, gain, 1);
01947 
01948         cpl_image_delete(skylocalmap); skylocalmap = NULL;
01949 
01950         if (images) {
01951             if (standard) {
01952                 cpl_table *ext_table  = NULL;
01953                 cpl_table *flux_table = NULL;
01954 
01955                 ext_table = dfs_load_table(frameset, "EXTINCT_TABLE", 1);
01956                 flux_table = dfs_load_table(frameset, "STD_FLUX_TABLE", 1);
01957 
01958                 photcal = mos_photometric_calibration(images[0], 
01959                                                       startwavelength,
01960                                                       dispersion, gain,
01961                                                       alltime, ext_table,
01962                                                       airmass, flux_table,
01963                                                       res_order);
01964 
01965                 cpl_table_delete(ext_table);
01966                 cpl_table_delete(flux_table);
01967 
01968                 if (photcal) {
01969 
01970                     if (qc) {
01971 
01972                         float *data;
01973                         char  *pipefile = NULL;
01974                         char   keyname[30];
01975 
01976                         qclist = dfs_load_header(frameset, science_tag, 0);
01977 
01978                         if (qclist == NULL)
01979                             fors_science_exit("Cannot reload scientific "
01980                                               "frame header");
01981 
01982                         fors_qc_start_group(qclist, "2.0", instrume);
01983 
01984 
01985                         /*
01986                          * QC1 group header
01987                          */
01988 
01989                         if (fors_qc_write_string("PRO.CATG", specphot_tag,
01990                                                  "Product category", instrume))
01991                             fors_science_exit("Cannot write product category "
01992                                               "to QC log file");
01993 
01994                         if (fors_qc_keyword_to_paf(qclist, "ESO DPR TYPE", NULL,
01995                                                    "DPR type", instrume))
01996                             fors_science_exit("Missing keyword DPR TYPE in "
01997                                               "scientific frame header");
01998 
01999                         if (fors_qc_keyword_to_paf(qclist, "ESO TPL ID", NULL,
02000                                                    "Template", instrume))
02001                             fors_science_exit("Missing keyword TPL ID in "
02002                                               "scientific frame header");
02003 
02004                         if (fors_qc_keyword_to_paf(qclist, 
02005                                                    "ESO INS GRIS1 NAME", NULL,
02006                                                    "Grism name", instrume))
02007                             fors_science_exit("Missing keyword INS GRIS1 NAME "
02008                                               "in scientific frame header");
02009 
02010                         if (fors_qc_keyword_to_paf(qclist, 
02011                                                    "ESO INS GRIS1 ID", NULL,
02012                                                    "Grism identifier", 
02013                                                    instrume))
02014                             fors_science_exit("Missing keyword INS GRIS1 ID "
02015                                               "in scientific frame header");
02016 
02017                         if (cpl_propertylist_has(qclist, "ESO INS FILT1 NAME"))
02018                             fors_qc_keyword_to_paf(qclist, 
02019                                                    "ESO INS FILT1 NAME", NULL,
02020                                                    "Filter name", instrume);
02021 
02022                         if (fors_qc_keyword_to_paf(qclist, 
02023                                                    "ESO INS COLL NAME", NULL,
02024                                                    "Collimator name", 
02025                                                    instrume))
02026                             fors_science_exit("Missing keyword INS COLL NAME "
02027                                               "in scientific frame header");
02028 
02029                         if (fors_qc_keyword_to_paf(qclist, 
02030                                                    "ESO DET CHIP1 ID", NULL,
02031                                                    "Chip identifier", 
02032                                                    instrume))
02033                             fors_science_exit("Missing keyword DET CHIP1 ID "
02034                                               "in scientific frame header");
02035 
02036                         if (fors_qc_keyword_to_paf(qclist, 
02037                                                    "ESO INS MOS10 WID",
02038                                                    "arcsec", "Slit width", 
02039                                                    instrume)) {
02040                             cpl_error_reset();
02041                             if (fors_qc_keyword_to_paf(qclist, 
02042                                                        "ESO INS SLIT WID",
02043                                                        "arcsec", "Slit width", 
02044                                                        instrume)) {
02045                                 if (mos) {
02046                                     fors_science_exit("Missing keyword "
02047                                                   "ESO INS MOS10 WID in "
02048                                                   "scientific frame header");
02049                                 }
02050                                 else {
02051                                     fors_science_exit("Missing keyword "
02052                                                   "ESO INS SLIT WID in "
02053                                                   "scientific frame header");
02054                                 }
02055                             }
02056                         }
02057 
02058 /*
02059                         if (mos) {
02060                             if (fors_qc_keyword_to_paf(qclist, 
02061                                                        "ESO INS MOS10 WID",
02062                                                        "arcsec", "Slit width", 
02063                                                        instrume))
02064                                 fors_science_exit("Missing keyword "
02065                                                   "ESO INS MOS10 WID in "
02066                                                   "scientific frame header");
02067                         }
02068                         else {
02069                             if (fors_qc_keyword_to_paf(qclist, 
02070                                                        "ESO INS SLIT WID",
02071                                                        "arcsec", "Slit width", 
02072                                                        instrume))
02073                                 fors_science_exit("Missing keyword "
02074                                                   "ESO INS SLIT WID in "
02075                                                   "scientific frame header");
02076                         }
02077 */
02078 
02079                         if (fors_qc_keyword_to_paf(qclist, 
02080                                                    "ESO DET WIN1 BINX", NULL,
02081                                                    "Binning factor along X", 
02082                                                    instrume))
02083                             fors_science_exit("Missing keyword ESO "
02084                                               "DET WIN1 BINX "
02085                                               "in scientific frame header");
02086         
02087                         if (fors_qc_keyword_to_paf(qclist, 
02088                                                    "ESO DET WIN1 BINY", NULL,
02089                                                    "Binning factor along Y", 
02090                                                    instrume))
02091                             fors_science_exit("Missing keyword "
02092                                               "ESO DET WIN1 BINY "
02093                                               "in scientific frame header");
02094 
02095                         if (fors_qc_keyword_to_paf(qclist, "ARCFILE", NULL,
02096                                                    "Archive name of input data",
02097                                                    instrume))
02098                             fors_science_exit("Missing keyword ARCFILE in "
02099                                               "scientific frame header");
02100         
02101                         pipefile = dfs_generate_filename_tfits(specphot_tag);
02102                         if (fors_qc_write_string("PIPEFILE", pipefile,
02103                                                  "Pipeline product name", 
02104                                                  instrume))
02105                             fors_science_exit("Cannot write PIPEFILE to "
02106                                               "QC log file");
02107                         cpl_free(pipefile);
02108 
02109 
02110                         /*
02111                          * QC1 parameters
02112                          */
02113 
02114                         wstart = 3700.;
02115                         wstep  = 400.;
02116                         wcount = 15;
02117 
02118                         dummy = cpl_image_new(wcount, 1, CPL_TYPE_FLOAT);
02119                         data = cpl_image_get_data_float(dummy);
02120                         map_table(dummy, wstart, wstep, photcal, 
02121                                   "WAVE", "EFFICIENCY");
02122 
02123                         for (i = 0; i < wcount; i++) {
02124                             sprintf(keyname, "QC.SPEC.EFFICIENCY%d.LAMBDA", 
02125                                     i + 1);
02126                             if (fors_qc_write_qc_double(qclist, 
02127                                                         wstart + wstep * i,
02128                                                         keyname, "Angstrom",
02129                                                         "Wavelength of "
02130                                                         "efficiency evaluation",
02131                                                         instrume)) {
02132                                 fors_science_exit("Cannot write wavelength of "
02133                                                   "efficiency evaluation");
02134                             }
02135 
02136                             sprintf(keyname, "QC.SPEC.EFFICIENCY%d", i + 1);
02137                             if (fors_qc_write_qc_double(qclist,
02138                                                         data[i],
02139                                                         keyname, "e-/photon",
02140                                                         "Efficiency",
02141                                                         instrume)) {
02142                                 fors_science_exit("Cannot write wavelength of "
02143                                                   "efficiency evaluation");
02144                             }
02145                         }
02146 
02147                         cpl_image_delete(dummy); dummy = NULL;
02148 
02149                         fors_qc_end_group();
02150         
02151                     }  /* End of QC1 computation */
02152 
02153                     if (dfs_save_table(frameset, photcal, specphot_tag, qclist,
02154                                        parlist, recipe, version))
02155                         fors_science_exit(NULL);
02156 
02157                     cpl_propertylist_delete(qclist); qclist = NULL;
02158 
02159                     if (have_phot) {
02160                         cpl_table_delete(photcal); photcal = NULL;
02161                     }
02162                 }
02163             }
02164 
02165             if (photometry) {
02166                 cpl_image *calibrated;
02167                 cpl_table *ext_table;
02168 
02169                 ext_table = dfs_load_table(frameset, "EXTINCT_TABLE", 1);
02170 
02171                 if (have_phot) {
02172                     if (cpl_frameset_count_tags(frameset, specphot_tag) == 0) {
02173                         photcal = dfs_load_table(frameset, 
02174                                                  master_specphot_tag, 1);
02175                     }
02176                     else {
02177                         photcal = dfs_load_table(frameset, specphot_tag, 1);
02178                     }
02179                 }
02180 
02181                 calibrated = mos_apply_photometry(images[0], photcal, 
02182                                                   ext_table, startwavelength, 
02183                                                   dispersion, gain, alltime, 
02184                                                   airmass);
02185                 cpl_propertylist_update_string(header, "BUNIT", 
02186                                    "10^(-16) erg/(cm^2 s Angstrom)");
02187 
02188                 if (dfs_save_image(frameset, calibrated,
02189                                    reduced_flux_science_tag, header,
02190                                    parlist, recipe, version)) {
02191                     cpl_image_delete(calibrated);
02192                     fors_science_exit(NULL);
02193                 }
02194 
02195                 cpl_table_delete(ext_table);
02196                 cpl_image_delete(calibrated);
02197             }
02198 
02199             if (time_normalise) {
02200                 cpl_propertylist_update_string(header, "BUNIT", "ADU/s");
02201                 cpl_image_divide_scalar(images[0], alltime);
02202             }
02203             else {
02204                 cpl_propertylist_update_string(header, "BUNIT", "ADU");
02205             }
02206 
02207             if (dfs_save_image(frameset, images[0], reduced_science_tag, header,
02208                                parlist, recipe, version))
02209                 fors_science_exit(NULL);
02210 
02211             if (time_normalise)
02212                 cpl_image_divide_scalar(images[1], alltime);
02213 
02214             if (dfs_save_image(frameset, images[1], reduced_sky_tag, header,
02215                                parlist, recipe, version))
02216                 fors_science_exit(NULL);
02217 
02218             cpl_image_delete(images[1]);
02219 
02220             if (photometry) {
02221                 cpl_image *calibrated;
02222                 cpl_table *ext_table; 
02223  
02224                 ext_table = dfs_load_table(frameset, "EXTINCT_TABLE", 1);
02225 
02226                 calibrated = mos_propagate_photometry_error(images[0], 
02227                                                   images[2], photcal,
02228                                                   ext_table, startwavelength,
02229                                                   dispersion, gain, alltime,
02230                                                   airmass);
02231 
02232                 cpl_propertylist_update_string(header, "BUNIT", 
02233                                    "10^(-16) erg/(cm^2 s Angstrom)");
02234 
02235                 if (dfs_save_image(frameset, calibrated,
02236                                    reduced_flux_error_tag, header,
02237                                    parlist, recipe, version)) {
02238                     cpl_image_delete(calibrated);
02239                     fors_science_exit(NULL);
02240                 }
02241 
02242                 cpl_table_delete(ext_table);
02243                 cpl_image_delete(calibrated);
02244             }
02245 
02246     
02247             if (time_normalise) {
02248                 cpl_propertylist_update_string(header, "BUNIT", "ADU/s");
02249                 cpl_image_divide_scalar(images[2], alltime);
02250             }
02251             else {
02252                 cpl_propertylist_update_string(header, "BUNIT", "ADU");
02253             }
02254 
02255             if (dfs_save_image(frameset, images[2], reduced_error_tag, header,
02256                                parlist, recipe, version))
02257                 fors_science_exit(NULL);
02258 
02259             cpl_image_delete(images[0]);
02260             cpl_image_delete(images[2]);
02261 
02262             cpl_free(images);
02263         }
02264         else {
02265             cpl_msg_warning(recipe, "No objects found: the products "
02266                             "%s, %s, and %s are not created", 
02267                             reduced_science_tag, reduced_sky_tag, 
02268                             reduced_error_tag);
02269         }
02270 
02271     }
02272 
02273     cpl_free(instrume); instrume = NULL;
02274     cpl_table_delete(slits); slits = NULL;
02275 
02276     if (skyalign >= 0) {
02277         if (dfs_save_table(frameset, idscoeff, disp_coeff_sky_tag, NULL, 
02278                            parlist, recipe, version))
02279             fors_science_exit(NULL);
02280     }
02281 
02282     cpl_table_delete(idscoeff); idscoeff = NULL;
02283 
02284 /*    if (skyglobal || skymedian >= 0) {   NSS */
02285     if (skyglobal || skymedian || skylocal) {
02286 
02287         if (photometry && photcal) {
02288             cpl_image *calibrated;
02289             cpl_table *ext_table; 
02290 
02291             ext_table = dfs_load_table(frameset, "EXTINCT_TABLE", 1);
02292 
02293             calibrated = mos_apply_photometry(mapped, photcal,
02294                                               ext_table, startwavelength,
02295                                               dispersion, gain, alltime,
02296                                               airmass);
02297 
02298             cpl_propertylist_update_string(header, "BUNIT", 
02299                                "10^(-16) erg/(cm^2 s Angstrom)");
02300 
02301             if (dfs_save_image(frameset, calibrated,
02302                                mapped_flux_science_tag, header,
02303                                parlist, recipe, version)) {
02304                 cpl_image_delete(calibrated);
02305                 fors_science_exit(NULL);
02306             }
02307 
02308             cpl_table_delete(ext_table);
02309             cpl_image_delete(calibrated);
02310         }
02311 
02312         if (time_normalise) {
02313             cpl_propertylist_update_string(header, "BUNIT", "ADU/s");
02314             cpl_image_divide_scalar(mapped, alltime);
02315         }
02316         else {
02317             cpl_propertylist_update_string(header, "BUNIT", "ADU");
02318         }
02319 
02320         if (dfs_save_image(frameset, mapped, mapped_science_tag, header, 
02321                            parlist, recipe, version))
02322             fors_science_exit(NULL);
02323     }
02324 
02325     cpl_table_delete(photcal); photcal = NULL;
02326     cpl_image_delete(mapped); mapped = NULL;
02327     cpl_propertylist_delete(header); header = NULL;
02328 
02329     if (cpl_error_get_code()) {
02330         cpl_msg_error(cpl_error_get_where(), cpl_error_get_message());
02331         fors_science_exit(NULL);
02332     }
02333     else 
02334         return 0;
02335 }