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