KMOS Pipeline Reference Manual  1.2.5
kmo_std_star.c
00001 /* $Id: kmo_std_star.c,v 1.79 2013/10/08 14:55:01 erw Exp $
00002  *
00003  * This file is part of the KMOS Pipeline
00004  * Copyright (C) 2002,2003 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., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00019  */
00020 
00021 /*
00022  * $Author: erw $
00023  * $Date: 2013/10/08 14:55:01 $
00024  * $Revision: 1.79 $
00025  * $Name: HEAD $
00026  */
00027 
00028 #ifdef HAVE_CONFIG_H
00029 #include <config.h>
00030 #endif
00031 
00032 /*-----------------------------------------------------------------------------
00033  *                              Includes
00034  *----------------------------------------------------------------------------*/
00035 
00036 #include <math.h>
00037 #include <string.h>
00038 
00039 #include <cpl.h>
00040 #include "kmclipm_math.h"
00041 
00042 #include "kmo_constants.h"
00043 #include "kmo_cpl_extensions.h"
00044 #include "kmo_utils.h"
00045 #include "kmo_functions.h"
00046 #include "kmo_priv_std_star.h"
00047 #include "kmo_priv_fit_profile.h"
00048 #include "kmo_priv_extract_spec.h"
00049 #include "kmo_priv_functions.h"
00050 #include "kmo_dfs.h"
00051 #include "kmo_error.h"
00052 #include "kmo_debug.h"
00053 #include "kmo_priv_reconstruct.h"
00054 
00055 /*-----------------------------------------------------------------------------
00056  *              Lines (all vacuum)
00057  *
00058 # Line lists for use in telluric transmission correction for KMOS, etc.
00059 # based on OBA standard stars.
00060 #
00061 # 30/01/2013   NMFS
00062 #
00063 #
00064 # - H lines of the Paschen and Brackett series (perhaps also Pfund series
00065 #   at the very red edge of K band) will be most prominent for late-O to
00066 #   A types.
00067 #
00068 # - HeI lines in absorption are mostly for O types (with some dependence
00069 #   on luminosity class).
00070 # - HeII lines will only be relevant in the earliest O types
00071 # - HeI and HeII lines may also appear in emission.
00072 #
00073 # The note "weak - irrelevant?" indicates lines that tend to be much
00074 # weaker, and would only be discernable in R > 5000 spectra with very
00075 # high S/N ratio.  They may cause asymmetric wings for neighbouring
00076 # stronger features depending on the star/spectral type.
00077 # They are included here for completeness, but can probably be ignored
00078 # in the context of KMOS telluric calibration.
00079 #
00080 # It is important, however, to include the stronger HeI and HeII features,
00081 # experience with SINFONI shows they are frequently there (esp. in H band).
00082 #
00083 #
00084 # N.B.
00085 #  The H line list in this file is complete within the Iz - K coverage
00086 #  of KMOS (excluding the highest Pa, Br, abd Pf transitions, which
00087 #  become very weak).
00088 #  The He line for >= 1.0um is fairly complete (strongest common lines
00089 #  are included).
00090 #  HOWEVER: the He line list at < 1.0um, relevant for Iz band, is missing.
00091 #
00092 #
00093 # Useful references:
00094 #  Wallace et al. 2000, ApJ, 535, 325
00095 #  Wallace & Hinkle 1997,
00096 #  Meyer et al. 1998,
00097 #  Hanson et al. 2005, ApJS, 161, 154
00098 #
00099 # In the future: planned XShooter stellar library (PI S. Trager) will
00100 # provide a cross-check over the full Iz - K band, as well as allow the
00101 # identification of potentially missing He features in the range 0.8-1um.
00102 
00103  *-----------------------------------------------------------------------------*/
00104 const int         nr_lines_h     = 10;
00105 const double      lines_center_h[]      = {
00106                                     1.7001,     // HeI          // triplet
00107                                     // 1.52616,    // Br-19        // (weak - irrelevant?)
00108                                     1.53429,    // Br-18
00109                                     1.54400,    // Br-17
00110                                     1.55576,    // Br-16
00111                                     1.57018,    // Br-15
00112                                     1.58817,    // Br-14
00113                                     1.61105,    // Br-13
00114                                     1.64084,    // Br-12
00115                                     1.68077,    // Br-11
00116                                     1.73634    // Br-10
00117                                     // 1.6918,     // HeII     // weak
00118                                     // 1.81754,    // Br-epsilon   // (in band, but low transmission)
00119                                     // 1.87524     // Pa-alpha     // (out of H-band? useful for HK?)
00120                                    };
00121 const double      lines_width_h[]  = {
00122                                     0.025,    // HeI
00123                                     // 0.015,    // Br-19
00124                                     0.003,    // Br-18
00125                                     0.015,    // Br-17
00126                                     0.015,    // Br-16
00127                                     0.015,    // Br-15
00128                                     0.025,    // Br-14
00129                                     0.015,    // Br-13
00130                                     0.025,    // Br-12
00131                                     0.025,    // Br-11
00132                                     0.05    // Br-10
00133                                     // 0.015,    // HeII
00134                                     // 0.015,    // Br-epsilon
00135                                     // 0.015     // Pa-alpha
00136                                        };
00137 const int         nr_lines_k     = 2;
00138 const double      lines_center_k[]      = {
00139                                     // 1.94470,    // Br-delta                  // out of K-band
00140                                     // 2.0581,     // HeI          // singlet   // very faint, non detectable
00141                                     2.1120,     // HeI          // triplet
00142                                     //2.1132,     // HeI          // singlet
00143                                     // 2.1494,     // HeI        // (weak - irrelevant?)
00144                                     // 2.1607,     // HeI        // triplet (weak - irrelevant?)
00145                                     // 2.1617,     // HeI        // singlet (weak - irrelevant?)
00146                                     // 2.1641,     // HeI        // triplet (weak - irrelevant?)
00147                                     2.16569    // Br-gamma
00148                                     // 2.1815,     // HeI        // (weak - irrelevant?)
00149                                     // 2.1840,     // HeI        // (weak - irrelevant?)
00150 // wo ?!?                                    2.1885,     // HeII
00151                                     // 2.43087,    // Pf-20      // (weak - irrelevant?)
00152                                     // 2.44851,    // Pf-19      // (weak - irrelevant?)
00153                                     // 2.46949,    // Pf-18      // (weak - irrelevant?)
00154                                     // 2.49475     // Pf-17      // (weak - irrelevant?)  // out of band
00155                                    };
00156 const double      lines_width_k[]   = {
00157                                     // 0.015,    // Br-delta     // (out of K-band? useful for HK?)
00158   //                                  0.008,     // HeI          // singlet
00159                                     0.01,     // HeI          // triplet
00160                                     //0.0015,     // HeI          // singlet
00161                                     // 0.003,     // HeI        // (weak - irrelevant?)
00162                                     // 0.003,     // HeI        // triplet (weak - irrelevant?)
00163                                     // 0.003,     // HeI        // singlet (weak - irrelevant?)
00164                                     // 0.015,     // HeI        // triplet (weak - irrelevant?)
00165                                     0.015    // Br-gamma
00166                                     // 0.003,     // HeI        // (weak - irrelevant?)
00167                                     // 0.003,     // HeI        // (weak - irrelevant?)
00168                                     // 0.015,     // HeII
00169                                     // 0.015,    // Pf-20      // (weak - irrelevant?)
00170                                     // 0.015,    // Pf-19      // (weak - irrelevant?)
00171                                     // 0.015,    // Pf-18      // (weak - irrelevant?)
00172                                     // 0.015     // Pf-17      // (weak - irrelevant?)
00173                                        };
00174 const int         nr_lines_hk    = 12;
00175 const double      lines_center_hk[]     = {
00176                         // H
00177                                     1.7001,     // HeI          // triplet
00178 
00179                                     1.53429,    // Br-18
00180                                     1.54400,    // Br-17
00181                                     1.55576,    // Br-16
00182                                     1.57018,    // Br-15
00183                                     1.58817,    // Br-14
00184                                     1.61105,    // Br-13
00185                                     1.64084,    // Br-12
00186                                     1.68077,    // Br-11
00187                                     1.73634,    // Br-10
00188                         // K
00189                                     2.1120,     // HeI          // triplet
00190                                     2.16569     // Br-gamma
00191                                    };
00192 const double      lines_width_hk[]  = {
00193                         // H
00194                                         0.025,    // HeI
00195                                         0.003,    // Br-18
00196                                         0.015,    // Br-17
00197                                         0.015,    // Br-16
00198                                         0.015,    // Br-15
00199                                         0.025,    // Br-14
00200                                         0.015,    // Br-13
00201                                         0.025,    // Br-12
00202                                         0.025,    // Br-11
00203                                         0.05,     // Br-10
00204                         // K
00205                                         0.015,     // HeI          // triplet
00206                                         0.015    // Br-gamma
00207                                        };
00208 const int         nr_lines_iz    = 12;
00209 const double      lines_center_iz[]  = {
00210                                         0.84386,    // Pa-18
00211                                         0.84679,    // Pa-17
00212                                         0.85031,    // Pa-16
00213                                         0.85460,    // Pa-15
00214                                         0.85990,    // Pa-14
00215                                         0.86657,    // Pa-13
00216                                         0.87511,    // Pa-12
00217                                         0.88635,    // Pa-11
00218                                         0.90156,    // Pa-10
00219                                         0.92297,    // Pa-9
00220                                         0.95467,    // Pa-epsilon
00221                                         1.00501     // Pa-delta
00222                                        };
00223 const double      lines_width_iz[]  = {
00224                                          0.0008,     // Pa-18
00225                                         0.003225,    // Pa-17
00226                                         0.0039,     // Pa-16
00227                                         0.0048,     // Pa-15
00228                                         0.006,     // Pa-14
00229                                         0.0076,     // Pa-13
00230                                         0.001,     // Pa-12
00231                                         0.013,     // Pa-11
00232                                          0.01,     // Pa-10
00233                                          0.013,     // Pa-9
00234                                          0.02,     // Pa-epsilon
00235                                          0.025       // Pa-delta
00236                                        };
00237 const int         nr_lines_yj    = 7;
00238 const double      lines_center_yj[]     = {
00239                                     // 1.00501,    // Pa-delta     // (out of band?)
00240                                     1.08331,    // HeI
00241                                     1.09160,    // HeI
00242                                     1.09389,    // Pa-gamma
00243 
00244                                     1.19723,    // HeI
00245 
00246                                     1.28191,    // Pa-beta
00247                                     1.27882,    // HeI
00248                                     // 1.28495,    // HeI   // faint
00249                                     1.29720     // HeI
00250                                    };
00251 const double      lines_width_yj[]  = {
00252                                     // 0.015,    // Pa-delta     // (out of band?)
00253                                     .01,//0.005,    // HeI
00254                                     .01,//0.002,    // HeI
00255                                     0.02,    // Pa-gamma
00256 
00257                                     0.003,    // HeI
00258 
00259                                     0.02,    // Pa-beta
00260                                     0.0025,    // HeI
00261                                     // 0.0007,    // HeI
00262                                     0.002     // HeI
00263                                        };
00264 
00265 /*-----------------------------------------------------------------------------
00266  *                          Functions prototypes
00267  *----------------------------------------------------------------------------*/
00268 
00269 static int kmo_std_star_create(cpl_plugin *);
00270 static int kmo_std_star_exec(cpl_plugin *);
00271 static int kmo_std_star_destroy(cpl_plugin *);
00272 static int kmo_std_star(cpl_parameterlist *, cpl_frameset *);
00273 
00274 /*-----------------------------------------------------------------------------
00275  *                          Static variables
00276  *----------------------------------------------------------------------------*/
00277 
00278 static char kmo_std_star_description[] =
00279 "This recipe creates a telluric calibration frame and a PSF frame. It must be\n"
00280 "called after the kmo_illumination-recipe.\n"
00281 "Since there won’t be enough standard stars to observe for all IFUs in one ex-\n"
00282 "posure, one has to do several exposures in a way that there is at least one\n"
00283 "standard star and one sky exposure in each IFU. A internal data organiser will\n"
00284 "analyse the provided exposures and select the appropriate frames as follows:\n"
00285 "1. For each IFU the first standard star in the list of provided exposures is\n"
00286 "   taken. All subsequent standard star exposures for this IFU will be ignored\n"
00287 "2. A corresponding sky exposure will be chosen which will be as close in time\n"
00288 "   to the standard star exposure as possible.\n"
00289 "3. For any IFUs not containing a standard star and a sky exposure an empty\n"
00290 "   frame will be returned.\n"
00291 "\n"
00292 "BASIC PARAMETERS:\n"
00293 "-----------------\n"
00294 "--startype\n"
00295 "If this parameter is specified, the stored star types of the observed obejcts \n"
00296 "in the FITS headers are overridden. This value applies to all objects exa-\n"
00297 "mined in the input frames. Examples would be “A3I”, “G3IV” or “K0I”. The first\n"
00298 "letter defines the star type, the second letter the spectral class and the last\n"
00299 "letters the luminosity class.\n"
00300 "\n"
00301 "--magnitude\n"
00302 "If this parameter is specified, the stored magnitudes in the FITS headers are \n"
00303 "overridden. For HK two magnitudes for each H and K have to be specified. All \n"
00304 "other gratings just use a single magnitude. If two values are provided, they \n"
00305 "have to be separated with a comma. \n"
00306 "\n"
00307 "--fmethod\n"
00308 "The type of function that should be fitted spatially to the collapsed image.\n"
00309 "This fit is used to create a mask to extract the spectrum of the object. Valid\n"
00310 "values are “gauss” and “moffat”.\n"
00311 "\n"
00312 "--imethod\n"
00313 "The interpolation method used for reconstruction. As default 'CS' is selected.\n"
00314 "Note that no error spectra will be generated for this interpolation method.\n"
00315 "Select a nearest neighbour method otherwise\n"
00316 "\n"
00317 "--range\n"
00318 "The spectral range [um] to combine when collapsing the reconstructed cubes.\n"
00319 "\n"
00320 "--save_cubes\n"
00321 "Set to TRUE if the intermediate reconstructed cubes (eventually divided by "
00322 "illumination correction) should be saved as well. Default is FALSE.\n"
00323 "\n"
00324 "ADVANCED PARAMETERS\n"
00325 "-------------------\n"
00326 "--flux\n"
00327 "Specify if flux conservation should be applied.\n"
00328 "\n"
00329 "--neighborhoodRange\n"
00330 "Defines the range to search for neighbors during reconstruction\n"
00331 "\n"
00332 "--b_samples\n"
00333 "The number of samples in spectral direction for the reconstructed cube.\n"
00334 "Ideally this number should be greater than 2048, the detector size.\n"
00335 "\n"
00336 "--b_start\n"
00337 "--b_end\n"
00338 "Used to define manually the start and end wavelength for the reconstructed\n"
00339 "cube. By default the internally defined values are used.\n"
00340 "\n"
00341 "--cmethod\n"
00342 "Following methods of frame combination are available:\n"
00343 "   * 'ksigma' (Default)\n"
00344 "   An iterative sigma clipping. For each position all pixels in the spectrum\n"
00345 "   are examined. If they deviate significantly, they will be rejected according\n"
00346 "   to the conditions:\n"
00347 "       val > mean + stdev * cpos_rej\n"
00348 "   and\n"
00349 "       val < mean - stdev * cneg_rej\n"
00350 "   where --cpos_rej, --cneg_rej and --citer are the corresponding configuration\n"
00351 "   parameters. In the first iteration median and percentile level are used.\n"
00352 "\n"
00353 "   * 'median'\n"
00354 "   At each pixel position the median is calculated.\n"
00355 "\n"
00356 "   * 'average'\n"
00357 "   At each pixel position the average is calculated.\n"
00358 "\n"
00359 "   * 'sum'\n"
00360 "   At each pixel position the sum is calculated.\n"
00361 "\n"
00362 "   * 'min_max'\n"
00363 "   The specified number of minimum and maximum pixel values will be rejected.\n"
00364 "   --cmax and --cmin apply to this method.\n"
00365 "\n"
00366 "--cpos_rej\n"
00367 "--cneg_rej\n"
00368 "--citer\n"
00369 "see --cmethod='ksigma'\n"
00370 "\n"
00371 "--cmax\n"
00372 "--cmin\n"
00373 "see --cmethod='min_max'\n"
00374 "\n"
00375 "--xcal_interpolation\n"
00376 "If true interpolate the pixel position in the slitlet (xcal) using the two\n"
00377 "closest rotator angles in the calibration file. Otherwise take the values\n"
00378 "of the closest rotator angle\n"
00379 "\n"
00380 "--suppress_extension\n"
00381 "If set to TRUE, the arbitrary filename extensions are supressed. If multiple\n"
00382 "products with the same category are produced, they will be numered consecutively\n"
00383 "starting from 0.\n"
00384 "\n"
00385 "-------------------------------------------------------------------------------\n"
00386 "  Input files:\n"
00387 "\n"
00388 "   DO                      KMOS                                                \n"
00389 "   category                Type  Explanation                   Required #Frames\n"
00390 "   --------                ----- -----------                   -------- -------\n"
00391 "   STD                     RAW   Std. star & sky exposures         Y     >=1   \n"
00392 "   XCAL                    F2D   x calibration frame               Y      1    \n"
00393 "   YCAL                    F2D   y calibration frame               Y      1    \n"
00394 "   LCAL                    F2D   Wavelength calib. frame           Y      1    \n"
00395 "   MASTER_FLAT             F2D   Master flat frame                 Y      1    \n"
00396 "   WAVE_BAND               F2L   Table with start-/end-wavelengths Y      1    \n"
00397 "   ILLUM_CORR              F2I   Illumination correction           N     0,1   \n"
00398 "   SOLAR_SPEC              F1S   Solar spectrum                    N     0,1   \n"
00399 "                                 (only for G stars)                            \n"
00400 "   ATMOS_MODEL             F1S   Model atmospheric transmisson     N     0,1   \n"
00401 "                                 (only for OBAF stars in K band)               \n"
00402 "   SPEC_TYPE_LOOKUP        F2L   LUT  eff. stellar temperature     N     0,1   \n"
00403 "\n"
00404 "  Output files:\n"
00405 "\n"
00406 "   DO                      KMOS\n"
00407 "   category                Type   Explanation\n"
00408 "   --------                -----  -----------\n"
00409 "   TELLURIC                F1I    The normalised telluric spectrum            \n"
00410 "                                  (including errors)                          \n"
00411 "   STAR_SPEC               F1I    The extracted star spectrum                 \n"
00412 "                                  (including errors)                          \n"
00413 "   STD_IMAGE               F2I    The standard star PSF images                \n"
00414 "   STD_MASK                F2I    The generated mask used to extract the star \n"
00415 "                                  spectrum                                    \n"
00416 "-------------------------------------------------------------------------------\n"
00417 "\n";
00418 
00419 /*-----------------------------------------------------------------------------
00420  *                              Functions code
00421  *----------------------------------------------------------------------------*/
00422 
00439 int cpl_plugin_get_info(cpl_pluginlist *list)
00440 {
00441     cpl_recipe *recipe = cpl_calloc(1, sizeof *recipe);
00442     cpl_plugin *plugin = &recipe->interface;
00443 
00444     cpl_plugin_init(plugin,
00445                         CPL_PLUGIN_API,
00446                         KMOS_BINARY_VERSION,
00447                         CPL_PLUGIN_TYPE_RECIPE,
00448                         "kmo_std_star",
00449                         "Create the telluric correction frame.",
00450                         kmo_std_star_description,
00451                         "Alex Agudo Berbel",
00452                         "kmos-spark@mpe.mpg.de",
00453                         kmos_get_license(),
00454                         kmo_std_star_create,
00455                         kmo_std_star_exec,
00456                         kmo_std_star_destroy);
00457 
00458     cpl_pluginlist_append(list, plugin);
00459 
00460     return 0;
00461 }
00462 
00470 static int kmo_std_star_create(cpl_plugin *plugin)
00471 {
00472     cpl_recipe *recipe;
00473     cpl_parameter *p;
00474 
00475     /* Check that the plugin is part of a valid recipe */
00476     if (cpl_plugin_get_type(plugin) == CPL_PLUGIN_TYPE_RECIPE) 
00477         recipe = (cpl_recipe *)plugin;
00478     else
00479         return -1;
00480 
00481     /* Create the parameters list in the cpl_recipe object */
00482     recipe->parameters = cpl_parameterlist_new();
00483 
00484     /* --startype */
00485     p = cpl_parameter_new_value("kmos.kmo_std_star.startype",
00486                                 CPL_TYPE_STRING,
00487                                 "The spectral type of the star (O, B, A, F, G)"
00488                                 " Format: G4V etc.",
00489                                 "kmos.kmo_std_star",
00490                                 "");
00491     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "startype");
00492     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00493     cpl_parameterlist_append(recipe->parameters, p);
00494 
00495     /* --imethod */
00496     p = cpl_parameter_new_value("kmos.kmo_std_star.imethod",
00497                                 CPL_TYPE_STRING,
00498                                 "Method to use for interpolation. "
00499                                 "[\"NN\" (nearest neighbour), "
00500                                 "\"lwNN\" (linear weighted nearest neighbor), "
00501                                 "\"swNN\" (square weighted nearest neighbor), "
00502                                 "\"MS\" (Modified Shepard's method), "
00503                                 "\"CS\" (Cubic spline)]",
00504                                 "kmos.kmo_std_star",
00505                                 "CS");
00506     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "imethod");
00507     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00508     cpl_parameterlist_append(recipe->parameters, p);
00509 
00510     /* --fmethod */
00511     p = cpl_parameter_new_value("kmos.kmo_std_star.fmethod",
00512                                 CPL_TYPE_STRING,
00513                                 "Either fit a 'gauss' or 'moffat' profile.",
00514                                 "kmos.kmo_std_star",
00515                                 "gauss");
00516     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "fmethod");
00517     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00518     cpl_parameterlist_append(recipe->parameters, p);
00519 
00520     /* --neighborhoodRange */
00521     p = cpl_parameter_new_value("kmos.kmo_std_star.neighborhoodRange",
00522                                 CPL_TYPE_DOUBLE,
00523                                 "Defines the range to search for neighbors "
00524                                 "in pixels",
00525                                 "kmos.kmo_std_star",
00526                                 1.001);
00527     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "neighborhoodRange");
00528     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00529     cpl_parameterlist_append(recipe->parameters, p);
00530 
00531     /* --magnitude */
00532     p = cpl_parameter_new_value("kmos.kmo_std_star.magnitude",
00533                                 CPL_TYPE_STRING,
00534                                 "The magnitude of the std star. For HK two "
00535                                 "values have to provided (eg. 12.1,13.2)",
00536                                 "kmos.kmo_std_star",
00537                                 "");
00538     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "magnitude");
00539     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00540     cpl_parameterlist_append(recipe->parameters, p);
00541 
00542     /* --flux */
00543     p = cpl_parameter_new_value("kmos.kmo_std_star.flux",
00544                                 CPL_TYPE_BOOL,
00545                                 "TRUE: Apply flux conservation. FALSE: otherwise",
00546                                 "kmos.kmo_std_star",
00547                                 TRUE);
00548     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "flux");
00549     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00550     cpl_parameterlist_append(recipe->parameters, p);
00551 
00552     /* --save_cubes */
00553     p = cpl_parameter_new_value("kmos.kmo_std_star.save_cubes",
00554                                 CPL_TYPE_BOOL,
00555                                 "TRUE: Save reconstructed cubes, FALSE: otherwise",
00556                                 "kmos.kmo_std_star",
00557                                 FALSE);
00558     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "save_cubes");
00559     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00560     cpl_parameterlist_append(recipe->parameters, p);
00561 
00562     /* --xcal_interpolation */
00563     p = cpl_parameter_new_value("kmos.kmo_std_star.xcal_interpolation",
00564                                 CPL_TYPE_BOOL,
00565                                 "TRUE: Interpolate xcal between rotator angles. FALSE: otherwise",
00566                                 "kmos.kmo_std_star",
00567                                 TRUE);
00568     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "xcal_interpolation");
00569     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00570     cpl_parameterlist_append(recipe->parameters, p);
00571 
00572     /* --suppress_extension */
00573     p = cpl_parameter_new_value("kmos.kmo_std_star.suppress_extension",
00574                                 CPL_TYPE_BOOL,
00575                                 "Suppress arbitrary filename extension."
00576                                 "(TRUE (apply) or FALSE (don't apply)",
00577                                 "kmos.kmo_std_star",
00578                                 FALSE);
00579     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "suppress_extension");
00580     cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
00581     cpl_parameterlist_append(recipe->parameters, p);
00582 
00583     // add parameters for band-definition
00584     kmo_band_pars_create(recipe->parameters,
00585                          "kmos.kmo_std_star");
00586 
00587     // add parameters for combining
00588     return kmo_combine_pars_create(recipe->parameters,
00589                                    "kmos.kmo_std_star",
00590                                    DEF_REJ_METHOD,
00591                                    FALSE);
00592 }
00593 
00599 static int kmo_std_star_exec(cpl_plugin *plugin)
00600 {
00601     cpl_recipe  *recipe;
00602 
00603     /* Get the recipe out of the plugin */
00604     if (cpl_plugin_get_type(plugin) == CPL_PLUGIN_TYPE_RECIPE) 
00605         recipe = (cpl_recipe *)plugin;
00606     else return -1;
00607 
00608     return kmo_std_star(recipe->parameters, recipe->frames);
00609 }
00610 
00616 static int kmo_std_star_destroy(cpl_plugin *plugin)
00617 {
00618     cpl_recipe *recipe;
00619 
00620     /* Get the recipe out of the plugin */
00621     if (cpl_plugin_get_type(plugin) == CPL_PLUGIN_TYPE_RECIPE) 
00622         recipe = (cpl_recipe *)plugin;
00623     else return -1 ;
00624 
00625     cpl_parameterlist_delete(recipe->parameters);
00626     return 0 ;
00627 }
00628 
00643 static int kmo_std_star(cpl_parameterlist *parlist, cpl_frameset *frameset)
00644 {
00645     cpl_imagelist    **stored_data_cube     = NULL,
00646                      **stored_noise_cube    = NULL;
00647     cpl_image        **stored_psf_data      = NULL,
00648                      *illum_corr            = NULL,
00649                      **stored_mask          = NULL,
00650                      *lcal                  = NULL;
00651     cpl_frame        *xcal_frame            = NULL,
00652                      *ycal_frame            = NULL,
00653                      *lcal_frame            = NULL,
00654                      *flat_frame            = NULL,
00655                      *illum_frame           = NULL,
00656                      *obj_frame             = NULL,
00657                      *sky_frame             = NULL,
00658                      *tmp_frame             = NULL;
00659     cpl_vector       *solar_spec             = NULL,
00660                      *atmos_model            = NULL,
00661                      **stored_telluric_data  = NULL,
00662                      **stored_telluric_noise = NULL,
00663                      **stored_starspec_data  = NULL,
00664                      **stored_starspec_noise = NULL,
00665                      *tmp_spec_data          = NULL,
00666                      *spec_qc     = NULL,
00667                      *tmp_spec_noise         = NULL,
00668                      *identified_slices      = NULL,
00669                      *tmp_vec                = NULL,
00670                      *lambda_x               = NULL;
00671     int              ret_val                = 0,
00672                      nr_devices             = 0,
00673                      nr_exp                 = 0,
00674                      j                      = 0,
00675                      *bounds                = NULL,
00676                      ifu_nr                 = 0,
00677                      citer                  = 0,
00678                      cmax                   = 0,
00679                      cmin                   = 0,
00680                      line_warning           = FALSE,
00681                      nr_std_stars           = 0,
00682                      print_warning_once     = TRUE,
00683                      flux                   = FALSE,
00684                      background             = FALSE,
00685                      band_method            = 0,
00686                      save_cubes             = FALSE,
00687                      has_magnitude          = TRUE,
00688                      xcal_interpolation     = FALSE,
00689                      suppress_extension     = FALSE,
00690                      nr_split_mag           = 0,
00691                      i                      = 0,
00692                      l                      = 0,
00693                      gx                     = 0,
00694                      gy                     = 0,
00695                      k                      = 0;
00696     const int        *punused_ifus          = NULL;
00697     objSkyStruct     *obj_sky_struct        = NULL;
00698     double           *stored_qc_throughput  = NULL,
00699                      star_temperature       = 0.0,
00700                      neighborhoodRange      = 1.001,
00701                      cpos_rej               = 0.0,
00702                      cneg_rej               = 0.0,
00703                      zeropoint              = -1.0,
00704                      throughput_mean        = -1.0,
00705                      throughput_sdv         = -1.0,
00706                      std_trace              = -1.0,
00707                      counts1                = 0.0,
00708                      counts2                = 0.0,
00709                      magnitude1             = 0.0,
00710                      magnitude2             = 0.0,
00711                      exptime                = 0.,
00712                      cdelt3                 = 0.,
00713                      mean_data              = 0.,
00714                      mean_noise             = 0.,
00715                      *ptmp_spec_noise       = NULL,
00716                      crpix1                 = 0.,
00717                      crval1                 = 0.,
00718                      cdelt1                 = 0.;
00719     const double     *ptmp_spec_data        = NULL;
00720     cpl_propertylist *main_header_tel                    = NULL,
00721                      *main_header_psf                    = NULL,
00722                      *sub_header_orig                    = NULL,
00723                      *tmp_sub_header                     = NULL,
00724                      *tmp_header                         = NULL,
00725                      **stored_sub_tel_data_headers       = NULL,
00726                      **stored_sub_tel_noise_headers      = NULL,
00727                      **stored_sub_cube_data_headers      = NULL,
00728                      **stored_sub_cube_noise_headers     = NULL,
00729                      **stored_sub_psf_headers            = NULL,
00730                      *pl_psf                             = NULL;
00731     cpl_table        *spec_type_LUT         = NULL,
00732                      *band_table            = NULL;;
00733     main_fits_desc   desc1,
00734                      desc2;
00735     char             *extname               = NULL,
00736                      *keyword               = NULL,
00737                      filename_telluric[256],
00738                      filename_starspec[256],
00739                      filename_psf[256],
00740                      filename_mask[256],
00741                      filename_cubes[256],
00742                      *suffix                = NULL,
00743                      *fn_suffix             = NULL,
00744                      spec_class[256],
00745                      lum_class[256],
00746                      star_type[2],
00747                      *tmp_band_method       = getenv("KMO_BAND_METHOD"),
00748                      **split_mag            = NULL,
00749                      *grat_id               = NULL;
00750     const char       *filter_id             = NULL,
00751                      *spec_type             = NULL,
00752                      *magnitude_txt         = NULL,
00753                      *imethod               = NULL,
00754                      *cmethod               = NULL,
00755                      *fmethod               = NULL,
00756                      *tmp_str               = NULL;
00757     gridDefinition   gd;
00758     cpl_array        **unused_ifus_before   = NULL,
00759                      **unused_ifus_after    = NULL;
00760     cpl_frameset     *frameset_std          = NULL;
00761 
00762     KMO_TRY
00763     {
00764         kmo_init_fits_desc(&desc1);
00765         kmo_init_fits_desc(&desc2);
00766 
00767         /* --- check input --- */
00768         KMO_TRY_ASSURE((parlist != NULL) &&
00769                        (frameset != NULL),
00770                        CPL_ERROR_NULL_INPUT,
00771                        "Not all input data is provided!");
00772 
00773         KMO_TRY_ASSURE(cpl_frameset_count_tags(frameset, STD) >= 1,
00774                        CPL_ERROR_ILLEGAL_INPUT,
00775                        "At least one STD frame is required!");
00776         if (cpl_frameset_count_tags(frameset, STD) == 1) {
00777             cpl_msg_warning("", "At least two STD frames should be provided "
00778                                 "in order to apply sky subtraction!");
00779         }
00780 
00781         KMO_TRY_ASSURE((cpl_frameset_count_tags(frameset, ILLUM_CORR) == 1) ||
00782                        (cpl_frameset_count_tags(frameset, ILLUM_CORR) == 0),
00783                        CPL_ERROR_FILE_NOT_FOUND,
00784                        "Exactly one ILLUM_CORR frame is required!");
00785         KMO_TRY_ASSURE(cpl_frameset_count_tags(frameset, XCAL) == 1,
00786                        CPL_ERROR_FILE_NOT_FOUND,
00787                        "Exactly one XCAL frame is required!");
00788         KMO_TRY_ASSURE(cpl_frameset_count_tags(frameset, YCAL) == 1,
00789                        CPL_ERROR_FILE_NOT_FOUND,
00790                        "Exactly one YCAL frame is required!");
00791         KMO_TRY_ASSURE(cpl_frameset_count_tags(frameset, LCAL) == 1,
00792                        CPL_ERROR_FILE_NOT_FOUND,
00793                        "Exactly one LCAL frame is required!");
00794         KMO_TRY_ASSURE(cpl_frameset_count_tags(frameset, MASTER_FLAT) == 1,
00795                        CPL_ERROR_FILE_NOT_FOUND,
00796                        "Exactly one MASTER_FLAT frame is required!");
00797         KMO_TRY_ASSURE(cpl_frameset_count_tags(frameset, WAVE_BAND) == 1,
00798                        CPL_ERROR_FILE_NOT_FOUND,
00799                        "Exactly one WAVE_BAND frame is required!");
00800         KMO_TRY_ASSURE(kmo_dfs_set_groups(frameset, "kmo_std_star") == 1,
00801                        CPL_ERROR_ILLEGAL_INPUT,
00802                        "Cannot identify RAW and CALIB frames!");
00803 
00804         /* --- get parameters --- */
00805         cpl_msg_info("", "--- Parameter setup for kmo_std_star ------");
00806 
00807         KMO_TRY_EXIT_IF_NULL(
00808             spec_type = kmo_dfs_get_parameter_string(parlist,
00809                                                      "kmos.kmo_std_star.startype"));
00810         KMO_TRY_EXIT_IF_ERROR(
00811             kmo_dfs_print_parameter_help(parlist,
00812                                          "kmos.kmo_std_star.startype"));
00813 
00814         KMO_TRY_EXIT_IF_NULL(
00815             imethod = kmo_dfs_get_parameter_string(parlist,
00816                                                    "kmos.kmo_std_star.imethod"));
00817         KMO_TRY_ASSURE((strcmp(imethod, "NN") == 0) ||
00818                        (strcmp(imethod, "lwNN") == 0) ||
00819                        (strcmp(imethod, "swNN") == 0) ||
00820                        (strcmp(imethod, "MS") == 0) ||
00821                        (strcmp(imethod, "CS") == 0),
00822                        CPL_ERROR_ILLEGAL_INPUT,
00823                        "method must be either \"NN\", \"lwNN\", "
00824                        "\"swNN\", \"MS\" or \"CS\"!");
00825         KMO_TRY_EXIT_IF_ERROR(
00826             kmo_dfs_print_parameter_help(parlist,
00827                                          "kmos.kmo_std_star.imethod"));
00828 
00829         KMO_TRY_EXIT_IF_NULL(
00830             fmethod = kmo_dfs_get_parameter_string(parlist,
00831                                                    "kmos.kmo_std_star.fmethod"));
00832         KMO_TRY_ASSURE((strcmp(fmethod, "gauss") == 0) ||
00833                        (strcmp(fmethod, "moffat") == 0),
00834                        CPL_ERROR_ILLEGAL_INPUT,
00835                        "fmethod must be either 'gauss' or "
00836                        "'moffat' !");
00837         KMO_TRY_EXIT_IF_ERROR(
00838             kmo_dfs_print_parameter_help(parlist,
00839                                         "kmos.kmo_std_star.method"));
00840 
00841         neighborhoodRange = kmo_dfs_get_parameter_double(parlist,
00842                                                          "kmos.kmo_std_star.neighborhoodRange");
00843         KMO_TRY_CHECK_ERROR_STATE();
00844         KMO_TRY_ASSURE(neighborhoodRange > 0.0,
00845                        CPL_ERROR_ILLEGAL_INPUT,
00846                        "neighborhoodRange must be greater than 0.0");
00847         KMO_TRY_EXIT_IF_ERROR(
00848             kmo_dfs_print_parameter_help(parlist,
00849                                          "kmos.kmo_std_star.neighborhoodRange"));
00850 
00851         magnitude_txt = kmo_dfs_get_parameter_string(parlist,
00852                                                      "kmos.kmo_std_star.magnitude");
00853         KMO_TRY_CHECK_ERROR_STATE();
00854         KMO_TRY_EXIT_IF_ERROR(
00855             kmo_dfs_print_parameter_help(parlist,
00856                                          "kmos.kmo_std_star.magnitude"));
00857 
00858         flux = kmo_dfs_get_parameter_bool(parlist,
00859                                           "kmos.kmo_std_star.flux");
00860         KMO_TRY_ASSURE((flux == FALSE) || (flux == TRUE),
00861                        CPL_ERROR_ILLEGAL_INPUT,
00862                        "flux must be either FALSE or TRUE!");
00863         KMO_TRY_EXIT_IF_ERROR(
00864             kmo_dfs_print_parameter_help(parlist, "kmos.kmo_std_star.flux"));
00865 
00866         save_cubes = kmo_dfs_get_parameter_bool(parlist,
00867                                           "kmos.kmo_std_star.save_cubes");
00868         KMO_TRY_ASSURE((save_cubes == FALSE) || (save_cubes == TRUE),
00869                        CPL_ERROR_ILLEGAL_INPUT,
00870                        "save_cubes must be either FALSE or TRUE!");
00871         KMO_TRY_EXIT_IF_ERROR(
00872             kmo_dfs_print_parameter_help(parlist, "kmos.kmo_std_star.save_cubes"));
00873 
00874         xcal_interpolation = kmo_dfs_get_parameter_bool(parlist,
00875                                            "kmos.kmo_std_star.xcal_interpolation");
00876         KMO_TRY_CHECK_ERROR_STATE();
00877         KMO_TRY_EXIT_IF_ERROR(
00878             kmo_dfs_print_parameter_help(parlist, "kmos.kmo_std_star.xcal_interpolation"));
00879         KMO_TRY_ASSURE((xcal_interpolation == TRUE) ||
00880                        (xcal_interpolation == FALSE),
00881                        CPL_ERROR_ILLEGAL_INPUT,
00882                        "xcal_interpolation must be TRUE or FALSE!");
00883 
00884         suppress_extension = kmo_dfs_get_parameter_bool(parlist,
00885                                           "kmos.kmo_std_star.suppress_extension");
00886         KMO_TRY_CHECK_ERROR_STATE();
00887         KMO_TRY_EXIT_IF_ERROR(
00888             kmo_dfs_print_parameter_help(parlist, "kmos.kmo_std_star.suppress_extension"));
00889 
00890         KMO_TRY_ASSURE((suppress_extension == TRUE) || (suppress_extension == FALSE),
00891                        CPL_ERROR_ILLEGAL_INPUT,
00892                        "suppress_extension must be TRUE or FALSE!");
00893 
00894         kmo_band_pars_load(parlist, "kmos.kmo_std_star");
00895 
00896         KMO_TRY_EXIT_IF_ERROR(
00897             kmo_combine_pars_load(parlist,
00898                                   "kmos.kmo_std_star",
00899                                   &cmethod,
00900                                   &cpos_rej,
00901                                   &cneg_rej,
00902                                   &citer,
00903                                   &cmin,
00904                                   &cmax,
00905                                   FALSE));
00906         cpl_msg_info("", "-------------------------------------------");
00907 
00908         //
00909         // Check if magnitude/frameset is valid and if throughput and zeropoint should be calculated
00910         //
00911 
00912         // Check if all STD frames have the same GRAT-ID
00913         // if not: don't calculate zeropoint and throughput
00914         KMO_TRY_EXIT_IF_NULL(
00915             frameset_std = cpl_frameset_new());
00916 
00917         KMO_TRY_EXIT_IF_NULL(
00918             tmp_frame = kmo_dfs_get_frame(frameset, STD));
00919         KMO_TRY_EXIT_IF_NULL(
00920             tmp_header = kmclipm_propertylist_load(cpl_frame_get_filename(tmp_frame), 0));
00921         KMO_TRY_EXIT_IF_NULL(
00922             grat_id = cpl_sprintf("%s", cpl_propertylist_get_string(tmp_header, "ESO INS GRAT1 ID")));
00923         KMO_TRY_EXIT_IF_ERROR(
00924             cpl_frameset_insert(frameset_std, cpl_frame_duplicate(tmp_frame)));
00925         cpl_propertylist_delete(tmp_header); tmp_header = NULL;
00926 KMO_TRY_CHECK_ERROR_STATE();
00927         KMO_TRY_EXIT_IF_NULL(
00928             tmp_frame = kmo_dfs_get_frame(frameset, NULL));
00929         while (tmp_frame != NULL ) {
00930             KMO_TRY_EXIT_IF_NULL(
00931                 tmp_header = kmclipm_propertylist_load(cpl_frame_get_filename(tmp_frame), 0));
00932             if (strcmp(grat_id, cpl_propertylist_get_string(tmp_header, "ESO INS GRAT1 ID")) == 0) {
00933                 // same grating
00934                 KMO_TRY_EXIT_IF_ERROR(
00935                     cpl_frameset_insert(frameset_std, cpl_frame_duplicate(tmp_frame)));
00936             } else {
00937                 // there are STD frames with different gratings
00938                 if (has_magnitude) {
00939                     cpl_msg_warning(cpl_func, "The STD frames have different gratings,"
00940                                             "following QC parameters won't be "
00941                                             "calculated: QC ZEROPOINT, QC THROUGHPUT,"
00942                                             "QC THROUGHPUT MEAN and QC THROUGHPUT STD");
00943                 }
00944                 has_magnitude = FALSE;
00945             }
00946             cpl_propertylist_delete(tmp_header); tmp_header = NULL;
00947 
00948             tmp_frame = kmo_dfs_get_frame(frameset, NULL);
00949             KMO_TRY_CHECK_ERROR_STATE();
00950         }
00951         KMO_TRY_CHECK_ERROR_STATE();
00952 
00953         if (cpl_frameset_count_tags(frameset, ATMOS_MODEL) == 1) {
00954             // check if ATMOS_MODEL is the band as the STD frames
00955             KMO_TRY_EXIT_IF_NULL(
00956                 tmp_frame = kmo_dfs_get_frame(frameset, ATMOS_MODEL));
00957             KMO_TRY_EXIT_IF_NULL(
00958                 tmp_sub_header = kmclipm_propertylist_load( cpl_frame_get_filename(tmp_frame), 0));
00959             KMO_TRY_EXIT_IF_NULL(
00960                 tmp_str = cpl_propertylist_get_string(tmp_sub_header, FILT_ID));
00961             KMO_TRY_ASSURE(strcmp(grat_id, tmp_str) == 0,
00962                            CPL_ERROR_ILLEGAL_INPUT,
00963                            "ATMOS model must have primary "
00964                            "keyword '%s' equal '%s'!!!",
00965                            FILT_ID, grat_id);
00966             cpl_propertylist_delete(tmp_sub_header);
00967             tmp_sub_header = NULL;
00968         }
00969 
00970         if (has_magnitude) {
00971             // all STD frames have the same GRAT-ID
00972             // now check source of magnitude (user or keyword)
00973             KMO_TRY_EXIT_IF_NULL(
00974                 tmp_frame = kmo_dfs_get_frame(frameset, STD));
00975             KMO_TRY_EXIT_IF_NULL(
00976                 tmp_header = kmclipm_propertylist_load(cpl_frame_get_filename(tmp_frame), 0));
00977 
00978             if (strcmp(magnitude_txt, "") == 0) {
00979                 // no user defined magnitude
00980 
00981                 // check for magnitude-keyword
00982                 if ((cpl_propertylist_has(tmp_header, STDSTAR_MAG)) &&
00983                     (cpl_propertylist_get_type(tmp_header, STDSTAR_MAG) == CPL_TYPE_STRING))
00984                 {
00985                     KMO_TRY_EXIT_IF_NULL(
00986                         magnitude_txt = cpl_propertylist_get_string(tmp_header, STDSTAR_MAG));
00987                     KMO_TRY_EXIT_IF_NULL(
00988                         split_mag = kmo_strsplit(magnitude_txt, ",", &nr_split_mag));
00989 
00990                     // check if band and number of magnitudes matches
00991                     if ((nr_split_mag == 2) &&
00992                         (strcmp(grat_id, "HK") == 0))
00993                     {
00994                         magnitude1 = atof(split_mag[0]);
00995                         magnitude2 = atof(split_mag[1]);
00996                         cpl_msg_info("", "Magnitude in H: %g", magnitude1);
00997                         cpl_msg_info("", "Magnitude in K: %g", magnitude2);
00998                     } else if ((nr_split_mag >= 1) &&
00999                                ((strcmp(grat_id, "K") == 0) ||
01000                                 (strcmp(grat_id, "H") == 0) ||
01001                                 (strcmp(grat_id, "IZ") == 0) ||
01002                                 (strcmp(grat_id, "YJ") == 0)))
01003                     {
01004                         magnitude1 = atof(split_mag[0]);
01005                         cpl_msg_info("", "Magnitude in %s: %g", grat_id, magnitude1);
01006                     } else {
01007                         // keyword STDSTAR_MAG doesn't match filter
01008                         has_magnitude = FALSE;
01009                         cpl_msg_warning(cpl_func, "The keyword %s doesn't match to grating',"
01010                                                   "following QC parameters won't be "
01011                                                   "calculated: QC ZEROPOINT, QC THROUGHPUT,"
01012                                                   "QC THROUGHPUT MEAN and QC THROUGHPUT STD", STDSTAR_MAG);
01013                     }
01014                     kmo_strfreev(split_mag);
01015                 } else {
01016                     // keyword STDSTAR_MAG unavailable or wrong type
01017                     has_magnitude = FALSE;
01018                     cpl_msg_warning(cpl_func, "The keyword %s is not set or of wrong type,"
01019                                               "following QC parameters won't be "
01020                                               "calculated: QC ZEROPOINT, QC THROUGHPUT,"
01021                                               "QC THROUGHPUT MEAN and QC THROUGHPUT STD", STDSTAR_MAG);
01022                 }
01023             } else {
01024                 // magnitude is user specified
01025                 cpl_msg_info(cpl_func, "Magnitude has been specified by user. Any "
01026                                        "value in keyword %s will be ignored.", STDSTAR_MAG);
01027 
01028                 KMO_TRY_EXIT_IF_NULL(
01029                     split_mag = kmo_strsplit(magnitude_txt, ",", &nr_split_mag));
01030                 switch (nr_split_mag) {
01031                 case 1:
01032                     magnitude1 = atof(split_mag[0]);
01033                     cpl_msg_info("", "Magnitude in %s: %g", grat_id, magnitude1);
01034                     break;
01035                 case 2:
01036                     magnitude1 = atof(split_mag[0]);
01037                     magnitude2 = atof(split_mag[1]);
01038                     cpl_msg_info("", "Magnitude in H: %g", magnitude1);
01039                     cpl_msg_info("", "Magnitude in K: %g", magnitude2);
01040                     break;
01041                 default:
01042                     KMO_TRY_ASSURE(1 == 0,
01043                                    CPL_ERROR_ILLEGAL_INPUT,
01044                                    "Provided magnitude was in wrong format! "
01045                                    "Either a single float value or two separated by comma");
01046                 }
01047                 kmo_strfreev(split_mag);
01048             }
01049             cpl_propertylist_delete(tmp_header); tmp_header = NULL;
01050         } // if (has_magnitude)
01051         cpl_msg_info("", "-------------------------------------------");
01052         KMO_TRY_CHECK_ERROR_STATE();
01053 
01054         //
01055         // check for spectral type (--startype) (user or keyword)
01056         //
01057         if (strcmp(spec_type, "") == 0) {
01058             // no user defined startype
01059 
01060             KMO_TRY_EXIT_IF_NULL(
01061                 tmp_frame = kmo_dfs_get_frame(frameset, STD));
01062             KMO_TRY_EXIT_IF_NULL(
01063                 tmp_header = kmclipm_propertylist_load(cpl_frame_get_filename(tmp_frame), 0));
01064 
01065             // check for startype-keyword
01066             if ((cpl_propertylist_has(tmp_header, STDSTAR_TYPE)) &&
01067                 (cpl_propertylist_get_type(tmp_header, STDSTAR_TYPE) == CPL_TYPE_STRING))
01068             {
01069                 KMO_TRY_EXIT_IF_NULL(
01070                     spec_type = cpl_propertylist_get_string(tmp_header, STDSTAR_TYPE));
01071             } else {
01072                 // keyword STDSTAR_TYPE unavailable or wrong type
01073             }
01074         } else {
01075             // startype is user specified
01076             cpl_msg_info(cpl_func, "Type of star has been specified by user. Any "
01077                                    "value in keyword %s will be ignored.", STDSTAR_TYPE);
01078         }
01079         KMO_TRY_CHECK_ERROR_STATE();
01080 
01081         if (strlen(spec_type) > 0) {
01082             if (kmo_get_spec_type(spec_type, spec_class, lum_class) != CPL_ERROR_NONE) {
01083                 cpl_error_reset();
01084                 spec_class[0] = '\0';
01085                 lum_class[0] = '\0';
01086                 star_type[0] = '\0';
01087                 cpl_msg_warning("", "The keyword %s is not set or of wrong type nor was it provided by the user. "
01088                                     "Can't divide solar spectrum for G stars or fit a profile "
01089                                     "to atmospheric transmission for OBAF stars and can't "
01090                                     "divide blackbody for any star.", STDSTAR_TYPE);
01091                 cpl_msg_warning("", "%s = '%s' (should be something like e.g.'G2V' odr 'A9III')", STDSTAR_TYPE, spec_type);
01092             } else {
01093                 strncpy(star_type, spec_class, 1);
01094                 star_type[1] = '\0';
01095                 cpl_msg_info("", "Spectral class:   %s", spec_class);
01096                 cpl_msg_info("", "Luminosity class: %s", lum_class);
01097             }
01098         } else {
01099             spec_class[0] = '\0';
01100             lum_class[0] = '\0';
01101             star_type[0] = '\0';
01102             cpl_msg_warning("", "The keyword %s is not set nor was it provided by the user. "
01103                                 "Can't divide solar spectrum for G stars or fit a profile "
01104                                 "to atmospheric transmission for OBAF stars and can't "
01105                                 "divide blackbody for any star.", STDSTAR_TYPE);
01106         }
01107         cpl_propertylist_delete(tmp_header); tmp_header = NULL;
01108         cpl_msg_info("", "-------------------------------------------");
01109         KMO_TRY_CHECK_ERROR_STATE();
01110 
01111         // assure that filters, grating and rotation offsets match for
01112         // XCAL, YCAL, LCAL and for data frame to reconstruct (except DARK
01113         // frames)
01114         // check if filter_id and grating_id match for all detectors
01115         KMO_TRY_EXIT_IF_ERROR(
01116             kmo_check_frameset_setup(frameset, XCAL, FALSE, FALSE, TRUE));
01117         KMO_TRY_EXIT_IF_ERROR(
01118             kmo_check_frame_setup(frameset, XCAL, YCAL, TRUE, FALSE, TRUE));
01119         KMO_TRY_EXIT_IF_ERROR(
01120             kmo_check_frame_setup(frameset, XCAL, LCAL, TRUE, FALSE, TRUE));
01121         KMO_TRY_EXIT_IF_ERROR(
01122             kmo_check_frame_setup(frameset, XCAL, MASTER_FLAT, TRUE, FALSE, TRUE));
01123         KMO_TRY_EXIT_IF_ERROR(
01124             kmo_check_frame_setup(frameset, XCAL, STD, FALSE, FALSE, TRUE));
01125 
01126         if (cpl_frameset_count_tags(frameset, ILLUM_CORR) == 1) {
01127             KMO_TRY_EXIT_IF_ERROR(
01128                 kmo_check_frame_setup(frameset, XCAL, ILLUM_CORR, TRUE, FALSE, FALSE));
01129         }
01130 
01131         // check descriptors of all frames
01132         KMO_TRY_EXIT_IF_NULL(
01133             xcal_frame = kmo_dfs_get_frame(frameset, XCAL));
01134 
01135         desc1 = kmo_identify_fits_header(cpl_frame_get_filename(xcal_frame));
01136         KMO_TRY_CHECK_ERROR_STATE();
01137 
01138         KMO_TRY_ASSURE((desc1.nr_ext % 3 == 0) &&
01139                        (desc1.ex_badpix == FALSE) &&
01140                        (desc1.fits_type == f2d_fits) &&
01141                        (desc1.frame_type == detector_frame),
01142                        CPL_ERROR_ILLEGAL_INPUT,
01143                        "XCAL isn't in the correct format!!!");
01144 
01145         KMO_TRY_EXIT_IF_NULL(
01146             ycal_frame = kmo_dfs_get_frame(frameset, YCAL));
01147         desc2 = kmo_identify_fits_header(cpl_frame_get_filename(ycal_frame));
01148         KMO_TRY_CHECK_ERROR_STATE();
01149 
01150         KMO_TRY_ASSURE((desc1.nr_ext == desc2.nr_ext) &&
01151                        (desc1.ex_badpix == desc2.ex_badpix) &&
01152                        (desc1.fits_type == desc2.fits_type) &&
01153                        (desc1.frame_type == desc2.frame_type),
01154                        CPL_ERROR_ILLEGAL_INPUT,
01155                        "YCAL isn't in the correct format!!!");
01156         kmo_free_fits_desc(&desc2);
01157         kmo_init_fits_desc(&desc2);
01158 
01159         KMO_TRY_EXIT_IF_NULL(
01160             lcal_frame = kmo_dfs_get_frame(frameset, LCAL));
01161         desc2 = kmo_identify_fits_header(cpl_frame_get_filename(lcal_frame));
01162         KMO_TRY_CHECK_ERROR_STATE();
01163 
01164         KMO_TRY_ASSURE((desc1.nr_ext == desc2.nr_ext) &&
01165                        (desc1.ex_badpix == desc2.ex_badpix) &&
01166                        (desc1.fits_type == desc2.fits_type) &&
01167                        (desc1.frame_type == desc2.frame_type),
01168                        CPL_ERROR_ILLEGAL_INPUT,
01169                        "YCAL isn't in the correct format!!!");
01170         kmo_free_fits_desc(&desc2);
01171         kmo_init_fits_desc(&desc2);
01172 
01173         KMO_TRY_EXIT_IF_NULL(
01174             flat_frame = kmo_dfs_get_frame(frameset, MASTER_FLAT));
01175         desc2 = kmo_identify_fits_header(cpl_frame_get_filename(flat_frame));
01176         KMO_TRY_CHECK_ERROR_STATE();
01177 
01178         KMO_TRY_ASSURE((desc2.nr_ext % 6 == 0) &&
01179                        (desc1.ex_badpix == desc2.ex_badpix) &&
01180                        (desc1.fits_type == desc2.fits_type) &&
01181                        (desc1.frame_type == desc2.frame_type),
01182                        CPL_ERROR_ILLEGAL_INPUT,
01183                        "MASTER_FLAT isn't in the correct format!!!");
01184         kmo_free_fits_desc(&desc2);
01185         kmo_init_fits_desc(&desc2);
01186 
01187         if (cpl_frameset_count_tags(frameset, ILLUM_CORR) == 1) {
01188             KMO_TRY_EXIT_IF_NULL(
01189                 illum_frame = kmo_dfs_get_frame(frameset, ILLUM_CORR));
01190             desc2 = kmo_identify_fits_header(cpl_frame_get_filename(illum_frame));
01191             KMO_TRY_CHECK_ERROR_STATE();
01192             KMO_TRY_ASSURE(((desc2.nr_ext == 24) || (desc2.nr_ext == 48)) &&
01193                            (desc2.ex_badpix == FALSE) &&
01194                            (desc2.fits_type == f2i_fits) &&
01195                            (desc2.frame_type == ifu_frame),
01196                            CPL_ERROR_ILLEGAL_INPUT,
01197                            "ILLUM_CORR isn't in the correct format!!!");
01198             kmo_free_fits_desc(&desc2);
01199             kmo_init_fits_desc(&desc2);
01200         }
01201 
01202         if (cpl_frameset_count_tags(frameset, SPEC_TYPE_LOOKUP) == 1) {
01203             KMO_TRY_EXIT_IF_NULL(
01204                 tmp_frame = kmo_dfs_get_frame(frameset, SPEC_TYPE_LOOKUP));
01205             desc2 = kmo_identify_fits_header(cpl_frame_get_filename(tmp_frame));
01206             KMO_TRY_CHECK_ERROR_STATE();
01207             KMO_TRY_ASSURE((desc2.nr_ext == 1) &&
01208                            (desc2.ex_badpix == FALSE) &&
01209                            (desc2.fits_type == f2l_fits) &&
01210                            (desc2.frame_type == list_frame),
01211                            CPL_ERROR_ILLEGAL_INPUT,
01212                            "SPEC_TYPE_LOOKUP isn't in the correct format!!!");
01213             kmo_free_fits_desc(&desc2);
01214             kmo_init_fits_desc(&desc2);
01215         }
01216 
01217         if (cpl_frameset_count_tags(frameset, SOLAR_SPEC) == 1) {
01218             KMO_TRY_EXIT_IF_NULL(
01219                 tmp_frame = kmo_dfs_get_frame(frameset, SOLAR_SPEC));
01220             desc2 = kmo_identify_fits_header(cpl_frame_get_filename(tmp_frame));
01221             KMO_TRY_CHECK_ERROR_STATE();
01222             KMO_TRY_ASSURE((desc2.nr_ext == 1) &&
01223                            (desc2.ex_badpix == FALSE) &&
01224                            (desc2.fits_type == f1s_fits) &&
01225                            (desc2.frame_type == spectrum_frame),
01226                            CPL_ERROR_ILLEGAL_INPUT,
01227                            "SOLAR_SPEC isn't in the correct format!!!");
01228             kmo_free_fits_desc(&desc2);
01229             kmo_init_fits_desc(&desc2);
01230         }
01231 
01232         if (cpl_frameset_count_tags(frameset, ATMOS_MODEL) == 1) {
01233             KMO_TRY_EXIT_IF_NULL(
01234                 tmp_frame = kmo_dfs_get_frame(frameset, ATMOS_MODEL));
01235             desc2 = kmo_identify_fits_header(cpl_frame_get_filename(tmp_frame));
01236             KMO_TRY_CHECK_ERROR_STATE();
01237             KMO_TRY_ASSURE((desc2.nr_ext == 1) &&
01238                            (desc2.ex_badpix == FALSE) &&
01239                            (desc2.fits_type == f1s_fits) &&
01240                            (desc2.frame_type == spectrum_frame),
01241                            CPL_ERROR_ILLEGAL_INPUT,
01242                            "ATMOS_MODEL isn't in the correct format!!!");
01243             kmo_free_fits_desc(&desc2);
01244             kmo_init_fits_desc(&desc2);
01245         }
01246 
01247         KMO_TRY_EXIT_IF_NULL(
01248             tmp_frame = kmo_dfs_get_frame(frameset, STD));
01249         while (tmp_frame != NULL ) {
01250             desc2 = kmo_identify_fits_header(cpl_frame_get_filename(tmp_frame));
01251             KMO_TRY_CHECK_ERROR_STATE();
01252             KMO_TRY_ASSURE((desc2.nr_ext == 3) &&
01253                            (desc2.ex_badpix == FALSE) &&
01254                            (desc2.fits_type == raw_fits) &&
01255                            (desc2.frame_type == detector_frame),
01256                            CPL_ERROR_ILLEGAL_INPUT,
01257                            "STD isn't in the correct format!!!");
01258             nr_devices = desc2.nr_ext;
01259             kmo_free_fits_desc(&desc2);
01260             kmo_init_fits_desc(&desc2);
01261 
01262             tmp_frame = kmo_dfs_get_frame(frameset, NULL);
01263             KMO_TRY_CHECK_ERROR_STATE();
01264         }
01265         KMO_TRY_EXIT_IF_NULL(
01266             tmp_frame = kmo_dfs_get_frame(frameset, STD));
01267         KMO_TRY_EXIT_IF_NULL(
01268             suffix = kmo_dfs_get_suffix(tmp_frame, TRUE, FALSE));
01269 
01270         KMO_TRY_EXIT_IF_ERROR(
01271             kmo_check_frame_setup_md5_xycal(frameset));
01272         KMO_TRY_EXIT_IF_ERROR(
01273             kmo_check_frame_setup_md5(frameset));
01274 
01275         cpl_msg_info("", "Detected instrument setup:   %s", suffix+1);
01276         cpl_msg_info("", "(grating 1, 2 & 3)");
01277 
01278         // check which IFUs are active for all frames
01279         KMO_TRY_EXIT_IF_NULL(
01280             unused_ifus_before = kmo_get_unused_ifus(frameset, 0, 0));
01281 
01282         KMO_TRY_EXIT_IF_NULL(
01283             unused_ifus_after = kmo_duplicate_unused_ifus(unused_ifus_before));
01284 
01285         kmo_print_unused_ifus(unused_ifus_before, FALSE);
01286 
01287         /* --- load data --- */
01288 
01289         if ((cpl_frameset_count_tags(frameset, SPEC_TYPE_LOOKUP) == 1) &&
01290             ((strlen(spec_class) > 0) || (strlen(lum_class) > 0)))
01291         {
01292             // get star temperature out of SPEC_TYPE_LOOKUP table
01293             KMO_TRY_EXIT_IF_NULL(
01294                 spec_type_LUT = kmo_dfs_load_table(frameset, SPEC_TYPE_LOOKUP, 1, 0));
01295             star_temperature = kmo_get_temperature(spec_type_LUT, spec_class, lum_class);
01296             KMO_TRY_CHECK_ERROR_STATE();
01297         } else if (cpl_frameset_count_tags(frameset, SPEC_TYPE_LOOKUP) != 1) {
01298             cpl_msg_warning("","No SPEC_TYPE_LOOKUP was provided! Can't divide blackbody.");
01299         } else if ((strlen(spec_class) == 0) || (strlen(lum_class) == 0)) {
01300 //            cpl_msg_warning("","No startype was provided! Can't "
01301 //                            "divide blackbody.");
01302         }
01303 
01304         // allocate intermediate memory
01305         KMO_TRY_EXIT_IF_NULL(
01306             stored_telluric_data    = (cpl_vector**)cpl_calloc(nr_devices*KMOS_IFUS_PER_DETECTOR,
01307                                                                sizeof(cpl_vector*)));
01308         KMO_TRY_EXIT_IF_NULL(
01309             stored_telluric_noise   = (cpl_vector**)cpl_calloc(nr_devices*KMOS_IFUS_PER_DETECTOR,
01310                                                                sizeof(cpl_vector*)));
01311         KMO_TRY_EXIT_IF_NULL(
01312             stored_starspec_data    = (cpl_vector**)cpl_calloc(nr_devices*KMOS_IFUS_PER_DETECTOR,
01313                                                                sizeof(cpl_vector*)));
01314         KMO_TRY_EXIT_IF_NULL(
01315             stored_starspec_noise   = (cpl_vector**)cpl_calloc(nr_devices*KMOS_IFUS_PER_DETECTOR,
01316                                                                sizeof(cpl_vector*)));
01317         KMO_TRY_EXIT_IF_NULL(
01318             stored_psf_data         = (cpl_image**)cpl_calloc(nr_devices*KMOS_IFUS_PER_DETECTOR,
01319                                                               sizeof(cpl_image*)));
01320         KMO_TRY_EXIT_IF_NULL(
01321             stored_mask             = (cpl_image**)cpl_calloc(nr_devices*KMOS_IFUS_PER_DETECTOR,
01322                                                               sizeof(cpl_image*)));
01323         KMO_TRY_EXIT_IF_NULL(
01324             stored_data_cube        = (cpl_imagelist**)cpl_calloc(nr_devices*KMOS_IFUS_PER_DETECTOR,
01325                                                               sizeof(cpl_imagelist*)));
01326         KMO_TRY_EXIT_IF_NULL(
01327             stored_noise_cube       = (cpl_imagelist**)cpl_calloc(nr_devices*KMOS_IFUS_PER_DETECTOR,
01328                                                               sizeof(cpl_imagelist*)));
01329         KMO_TRY_EXIT_IF_NULL(
01330             stored_qc_throughput    = (double*)cpl_calloc(nr_devices*KMOS_IFUS_PER_DETECTOR,
01331                                                           sizeof(double)));
01332         KMO_TRY_EXIT_IF_NULL(
01333             stored_sub_psf_headers  = (cpl_propertylist**)cpl_calloc(nr_devices*KMOS_IFUS_PER_DETECTOR,
01334                                                                      sizeof(cpl_propertylist*)));
01335         KMO_TRY_EXIT_IF_NULL(
01336             stored_sub_tel_data_headers = (cpl_propertylist**)cpl_calloc(nr_devices*KMOS_IFUS_PER_DETECTOR,
01337                                                                          sizeof(cpl_propertylist*)));
01338         KMO_TRY_EXIT_IF_NULL(
01339             stored_sub_tel_noise_headers = (cpl_propertylist**)cpl_calloc(nr_devices*KMOS_IFUS_PER_DETECTOR,
01340                                                                          sizeof(cpl_propertylist*)));
01341 
01342         if (save_cubes) {
01343             KMO_TRY_EXIT_IF_NULL(
01344                 stored_sub_cube_data_headers = (cpl_propertylist**)cpl_calloc(nr_devices*KMOS_IFUS_PER_DETECTOR,
01345                                                                              sizeof(cpl_propertylist*)));
01346             KMO_TRY_EXIT_IF_NULL(
01347                 stored_sub_cube_noise_headers = (cpl_propertylist**)cpl_calloc(nr_devices*KMOS_IFUS_PER_DETECTOR,
01348                                                                              sizeof(cpl_propertylist*)));
01349         }
01350 
01351         // get bounds
01352         KMO_TRY_EXIT_IF_NULL(
01353             tmp_header = kmo_dfs_load_primary_header(frameset, XCAL));
01354         KMO_TRY_EXIT_IF_NULL(
01355             bounds = kmclipm_extract_bounds(tmp_header));
01356         cpl_propertylist_delete(tmp_header); tmp_header = NULL;
01357 
01358         // setup grid definition, wavelength start and end points will be set
01359         // in the detector loop
01360         KMO_TRY_EXIT_IF_ERROR(
01361             kmclipm_setup_grid(&gd, imethod, neighborhoodRange, KMOS_PIX_RESOLUTION, 0.));
01362 
01363         // get valid STD frames with objects in it and associated sky exposures
01364         KMO_TRY_EXIT_IF_NULL(
01365             obj_sky_struct = kmo_create_objSkyStruct(frameset_std, STD, FALSE));
01366         kmo_print_objSkyStruct(obj_sky_struct);
01367 
01368         // loop the object-sky pairs
01369         if (obj_sky_struct->size == 0) {
01370             cpl_msg_warning(cpl_func,"Not a single frame contains an object");
01371         } else {
01372             strcpy(filename_telluric, TELLURIC);
01373             strcpy(filename_starspec, STAR_SPEC);
01374             strcpy(filename_psf, STD_IMAGE);
01375             strcpy(filename_mask, STD_MASK);
01376             strcpy(filename_cubes, STD_CUBE);
01377 
01378             obj_frame = obj_sky_struct->table[nr_exp].objFrame;
01379             KMO_TRY_EXIT_IF_NULL(
01380                 main_header_tel = kmclipm_propertylist_load(cpl_frame_get_filename(obj_frame), 0));
01381 
01382             exptime = cpl_propertylist_get_double(main_header_tel, EXPTIME);
01383             KMO_TRY_CHECK_ERROR_STATE();
01384 
01385             // load, process & store frames
01386 
01387             for (i = 1; i <= nr_devices; i++) {
01388                 // extract LCAL image close to ROTANGLE 0. assuming that the wavelength range
01389                 // doesn't differ too much with different ROTANGLEs.
01390                 print_cal_angle_msg_once = FALSE;
01391                 double rotangle_found;
01392                 KMO_TRY_EXIT_IF_NULL(
01393                     lcal = kmo_dfs_load_cal_image(frameset, LCAL, i, FALSE, 0., FALSE, NULL,
01394                             &rotangle_found, -1, 0, 0));
01395                 if (i==1) {
01396                     print_cal_angle_msg_once = TRUE;
01397                 }
01398                 if (tmp_band_method != NULL) {
01399                     band_method = atoi(tmp_band_method);
01400                 }
01401 
01402                 // get filter for this detector
01403                 // ESO INS FILTi ID
01404                 KMO_TRY_EXIT_IF_NULL(
01405                     keyword = cpl_sprintf("%s%d%s", IFU_FILTID_PREFIX, i, IFU_FILTID_POSTFIX));
01406                 filter_id = cpl_propertylist_get_string(main_header_tel, keyword);
01407                 cpl_free(keyword); keyword = NULL;
01408 
01409                 KMO_TRY_EXIT_IF_NULL(
01410                     band_table = kmo_dfs_load_table(frameset, WAVE_BAND, 1, 0));
01411                 KMO_TRY_EXIT_IF_ERROR(
01412                     kmclipm_setup_grid_band_lcal(&gd, lcal, filter_id,
01413                                                  band_method, band_table));
01414                 cpl_image_delete(lcal); lcal = NULL;
01415                 cpl_table_delete(band_table); band_table = NULL;
01416 
01417                 // load sub_header of original F2D image
01418                 KMO_TRY_EXIT_IF_NULL(
01419                     sub_header_orig = kmclipm_propertylist_load( cpl_frame_get_filename(obj_frame), i));
01420 
01421                 for (j = 0; j < KMOS_IFUS_PER_DETECTOR; j++) {
01422                     ifu_nr = (i-1)*KMOS_IFUS_PER_DETECTOR + j + 1;
01423                     // check if IFU is valid according to main header keywords &
01424                     // calibration files
01425                     // AND check if there is a sky frame available for this IFU
01426                     kmo_collapse_objSkyStruct(obj_sky_struct, ifu_nr,
01427                                               &obj_frame, &sky_frame);
01428 
01429                     KMO_TRY_EXIT_IF_NULL(
01430                         punused_ifus = cpl_array_get_data_int_const(unused_ifus_after[i-1]));
01431 
01432                     // Search for keyword ESO OCS ARMi NOTUSED
01433                     // If not present (CPL_ERROR_DATA_NOT_FOUND) we will eventually
01434                     // process standard star
01435                     KMO_TRY_EXIT_IF_NULL(
01436                         keyword = cpl_sprintf("%s%d%s", IFU_VALID_PREFIX, ifu_nr, IFU_VALID_POSTFIX));
01437                     tmp_str = cpl_propertylist_get_string(main_header_tel, keyword);
01438                     cpl_free(keyword); keyword = NULL;
01439 
01440                     if ((cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) &&
01441                         (bounds[2*(ifu_nr-1)] != -1) &&
01442                         (bounds[2*(ifu_nr-1)+1] != -1) &&
01443                         (sky_frame != NULL) &&
01444                         (punused_ifus[j] == 0))
01445                     {
01446                         cpl_error_reset();
01447                         // IFU is valid
01448 
01449                         if (sky_frame != NO_CORRESPONDING_SKYFRAME) {
01450                             cpl_msg_info("","Processing standard star in IFU %d "
01451                                          "(obj: %s, sky: %s)", ifu_nr,
01452                                          cpl_frame_get_filename(obj_frame),
01453                                          cpl_frame_get_filename(sky_frame));
01454                         } else {
01455                             sky_frame = NULL;
01456                             cpl_msg_warning("","Processing standard star in IFU %d "
01457                                          "(obj: %s, no corresponding sky frame",
01458                                          ifu_nr, cpl_frame_get_filename(obj_frame));
01459                         }
01460 
01461                         nr_std_stars++;
01462 
01463                         char *ggg = cpl_sprintf("%s%d", PRO_STD, ifu_nr);
01464                         KMO_TRY_EXIT_IF_ERROR(
01465                             cpl_propertylist_update_int(main_header_tel, ggg, 1));
01466                         cpl_free(ggg); ggg = NULL;
01467 
01468                         // calculate WCS and make copies of sub_header
01469                         KMO_TRY_EXIT_IF_NULL(
01470                             tmp_sub_header = cpl_propertylist_duplicate(sub_header_orig));
01471                         KMO_TRY_EXIT_IF_ERROR(
01472                             kmo_calc_wcs_gd(main_header_tel, tmp_sub_header, ifu_nr, gd));
01473                         KMO_TRY_EXIT_IF_NULL(
01474                             stored_sub_tel_data_headers[ifu_nr-1] =
01475                                    cpl_propertylist_duplicate(tmp_sub_header));
01476                         KMO_TRY_EXIT_IF_NULL(
01477                             stored_sub_psf_headers[ifu_nr-1] =
01478                                    cpl_propertylist_duplicate(tmp_sub_header));
01479                         if (save_cubes) {
01480                             KMO_TRY_EXIT_IF_NULL(
01481                                 stored_sub_cube_data_headers[ifu_nr-1] =
01482                                        cpl_propertylist_duplicate(tmp_sub_header));
01483                         }
01484                         cpl_propertylist_delete(tmp_sub_header);
01485                         tmp_sub_header = NULL;
01486 
01487                         //
01488                         // adjust telluric-headers: copy CRPIX3 to CRPIX1,
01489                         //
01490                         cpl_propertylist_update_double(stored_sub_tel_data_headers[ifu_nr-1], CRVAL1,
01491                                 cpl_propertylist_get_double(stored_sub_tel_data_headers[ifu_nr-1], CRVAL3));
01492                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CRVAL2);
01493                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CRVAL3);
01494                         KMO_TRY_CHECK_ERROR_STATE();
01495 
01496                         // CRPIX
01497                         cpl_propertylist_update_double(stored_sub_tel_data_headers[ifu_nr-1], CRPIX1,
01498                                 cpl_propertylist_get_double(stored_sub_tel_data_headers[ifu_nr-1], CRPIX3));
01499                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CRPIX2);
01500                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CRPIX3);
01501                         KMO_TRY_CHECK_ERROR_STATE();
01502 
01503                         // CDELT
01504                         cdelt3 = cpl_propertylist_get_double(stored_sub_tel_data_headers[ifu_nr-1], CDELT3);
01505                         cpl_propertylist_update_double(stored_sub_tel_data_headers[ifu_nr-1], CDELT1,
01506                                 cdelt3);
01507                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CDELT2);
01508                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CDELT3);
01509                         KMO_TRY_CHECK_ERROR_STATE();
01510 
01511                         // CTYPE
01512                         cpl_propertylist_update_string(stored_sub_tel_data_headers[ifu_nr-1], CTYPE1,
01513                                 cpl_propertylist_get_string(stored_sub_tel_data_headers[ifu_nr-1], CTYPE3));
01514                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CTYPE2);
01515                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CTYPE3);
01516                         KMO_TRY_CHECK_ERROR_STATE();
01517 
01518                         // CUNIT
01519                         cpl_propertylist_update_string(stored_sub_tel_data_headers[ifu_nr-1], CUNIT1,
01520                                 cpl_propertylist_get_string(stored_sub_tel_data_headers[ifu_nr-1], CUNIT3));
01521                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CUNIT2);
01522                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CUNIT3);
01523 
01524                         // CDx_x
01525                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CD1_1);
01526                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CD1_2);
01527                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CD1_3);
01528                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CD2_1);
01529                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CD2_2);
01530                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CD2_3);
01531                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CD3_1);
01532                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CD3_2);
01533                         cpl_propertylist_erase(stored_sub_tel_data_headers[ifu_nr-1], CD3_3);
01534                         KMO_TRY_CHECK_ERROR_STATE();
01535 
01536                         //
01537                         // adjust psf-headers: delete CRPIX3 etc.
01538                         //
01539                         cpl_propertylist_erase(stored_sub_psf_headers[ifu_nr-1], CRPIX3);
01540                         cpl_propertylist_erase(stored_sub_psf_headers[ifu_nr-1], CRPIX3);
01541                         cpl_propertylist_erase(stored_sub_psf_headers[ifu_nr-1], CDELT3);
01542                         cpl_propertylist_erase(stored_sub_psf_headers[ifu_nr-1], CRVAL3);
01543                         cpl_propertylist_erase(stored_sub_psf_headers[ifu_nr-1], CTYPE3);
01544                         cpl_propertylist_erase(stored_sub_psf_headers[ifu_nr-1], CUNIT3);
01545                         cpl_propertylist_erase(stored_sub_psf_headers[ifu_nr-1], CD1_3);
01546                         cpl_propertylist_erase(stored_sub_psf_headers[ifu_nr-1], CD2_3);
01547                         cpl_propertylist_erase(stored_sub_psf_headers[ifu_nr-1], CD3_1);
01548                         cpl_propertylist_erase(stored_sub_psf_headers[ifu_nr-1], CD3_2);
01549                         cpl_propertylist_erase(stored_sub_psf_headers[ifu_nr-1], CD3_3);
01550                         KMO_TRY_CHECK_ERROR_STATE();
01551 
01552                         KMO_TRY_EXIT_IF_ERROR(
01553                             kmo_reconstruct_sci(ifu_nr,
01554                                                 bounds[2*(ifu_nr-1)],
01555                                                 bounds[2*(ifu_nr-1)+1],
01556                                                 obj_frame,
01557                                                 STD,
01558                                                 sky_frame,
01559                                                 STD,
01560                                                 flat_frame,
01561                                                 xcal_frame,
01562                                                 ycal_frame,
01563                                                 lcal_frame,
01564                                                 NULL,
01565                                                 &gd,
01566                                                 &stored_data_cube[ifu_nr-1],
01567                                                 &stored_noise_cube[ifu_nr-1],
01568                                                 flux,
01569                                                 background,
01570                                                 xcal_interpolation));
01571 
01572                         // divide illumination correction from the data_cube
01573                         // (illumination noise will be very small versus
01574                         // noise_cube, so it is skipped here)
01575                         if (cpl_frameset_count_tags(frameset, ILLUM_CORR) == 1) {
01576                             KMO_TRY_EXIT_IF_NULL(
01577                                 illum_corr = kmo_dfs_load_image_frame(illum_frame, ifu_nr,
01578                                                                       FALSE, FALSE, NULL));
01579                             KMO_TRY_EXIT_IF_ERROR(
01580                                 cpl_imagelist_divide_image(stored_data_cube[ifu_nr-1], illum_corr));
01581                             cpl_image_delete(illum_corr); illum_corr = NULL;
01582                         }
01583 
01584                         // calculate QC_STD_TRACE
01585                         // (the distance of the PSF to the centre)
01586                         KMO_TRY_EXIT_IF_ERROR(
01587                             kmo_calculate_std_trace(stored_data_cube[ifu_nr-1], fmethod, &std_trace));
01588 
01589                         KMO_TRY_EXIT_IF_ERROR(
01590                             kmclipm_update_property_double(stored_sub_psf_headers[ifu_nr-1],
01591                                                            QC_STD_TRACE, std_trace,
01592                                                            "[pix] distance of PSF and centre of IFU"));
01593 
01594                         KMO_TRY_EXIT_IF_NULL(
01595                             identified_slices = cpl_vector_new(cpl_imagelist_get_size(stored_data_cube[ifu_nr-1])));
01596                         KMO_TRY_EXIT_IF_ERROR(
01597                             cpl_vector_fill(identified_slices, 1.0));
01598 
01599                         // collapse cube and get PSF image
01600                         KMO_TRY_EXIT_IF_ERROR(
01601                             kmclipm_make_image(stored_data_cube[ifu_nr-1], NULL,
01602                                                &stored_psf_data[ifu_nr-1], NULL,
01603                                                identified_slices,
01604                                                cmethod,
01605                                                cpos_rej, cneg_rej, citer,
01606                                                cmax, cmin));
01607                         cpl_vector_delete(identified_slices);
01608                         identified_slices= NULL;
01609 
01610                         // fit a 2D profile to get a mask and fwhm in x and y,
01611                         KMO_TRY_EXIT_IF_NULL(
01612                             tmp_vec = kmo_fit_profile_2D(stored_psf_data[ifu_nr-1],
01613                                                          NULL,
01614                                                          fmethod,
01615                                                          &stored_mask[ifu_nr-1],
01616                                                          &pl_psf));
01617 
01618                         // normalise mask to 1 and clip values below 0.5
01619                         cpl_image_divide_scalar(stored_mask[ifu_nr-1], cpl_image_get_max(stored_mask[ifu_nr-1]));
01620                         KMO_TRY_CHECK_ERROR_STATE();
01621 
01622                         int dummy=0;
01623                         for (gx = 1; gx <= cpl_image_get_size_x(stored_mask[ifu_nr-1]); gx++) {
01624                             for (gy = 1; gy <= cpl_image_get_size_y(stored_mask[ifu_nr-1]); gy++) {
01625                                 if (cpl_image_get(stored_mask[ifu_nr-1], gx, gy, &dummy) < 0.5) {
01626                                     cpl_image_set(stored_mask[ifu_nr-1], gx, gy, 0.);
01627                                 } else {
01628                                     cpl_image_set(stored_mask[ifu_nr-1], gx, gy, 1.);
01629                                 }
01630                             }
01631                         }
01632                         KMO_TRY_CHECK_ERROR_STATE();
01633 
01634                         // update subheader with fit parameters
01635                         KMO_TRY_EXIT_IF_ERROR(
01636                             cpl_propertylist_append(stored_sub_tel_data_headers[ifu_nr-1], pl_psf));
01637                         cpl_propertylist_delete(pl_psf); pl_psf = NULL;
01638 
01639                         // store QC_SPAT_RES (RMS of fwhm_x and fwhm_y)
01640                         double factor_fwhm = 2*sqrt(2*log(2));
01641                         double spat_res = pow(cpl_vector_get(tmp_vec, 4) * factor_fwhm, 2);
01642                         spat_res += pow(cpl_vector_get(tmp_vec, 5) * factor_fwhm, 2);
01643                         spat_res /= 2;
01644                         KMO_TRY_EXIT_IF_ERROR(
01645                             kmclipm_update_property_double(stored_sub_psf_headers[ifu_nr-1],
01646                                                            QC_SPAT_RES,
01647                                                            sqrt(spat_res)*KMOS_PIX_RESOLUTION,
01648                                                            "[arcsec] mean fwhm resolution of PSF"));
01649                         cpl_vector_delete(tmp_vec); tmp_vec = NULL;
01650 
01651                         // extract spectrum in masked area
01652                         KMO_TRY_EXIT_IF_ERROR(
01653                             kmo_priv_extract_spec(stored_data_cube[ifu_nr-1],
01654                                                   stored_noise_cube[ifu_nr-1],
01655                                                   stored_mask[ifu_nr-1],
01656                                                   &tmp_spec_data,
01657                                                   &tmp_spec_noise));
01658 
01659                         // store to save to disk later on
01660                         stored_starspec_data[ifu_nr-1] = cpl_vector_duplicate(tmp_spec_data);
01661                         if (tmp_spec_noise != NULL) {
01662                             stored_starspec_noise[ifu_nr-1] = cpl_vector_duplicate(tmp_spec_noise);
01663                         }
01664                         KMO_TRY_CHECK_ERROR_STATE();
01665 
01666                         // if magnitude is provided
01667                         // calculate zeropoint and throughput
01668                         if (has_magnitude) {
01669                             // extract spectrum of whole are for QC THRUHput and ZEROPOINT
01670                             KMO_TRY_EXIT_IF_ERROR(
01671                                 kmo_priv_extract_spec(stored_data_cube[ifu_nr-1],
01672                                                       NULL,
01673                                                       NULL,
01674                                                       &spec_qc,
01675                                                       NULL));
01676 
01677                             // multiply spectrum with area of IFU (196) to get the sum
01678                             const cpl_image *tmp_img = cpl_imagelist_get(stored_data_cube[ifu_nr-1], 0);
01679                             int tmpx = cpl_image_get_size_x(tmp_img),
01680                                 tmpy = cpl_image_get_size_x(tmp_img);
01681                             cpl_vector_multiply_scalar(spec_qc, tmpx*tmpy);
01682                         }
01683 
01684                         // calculate abscissa of output spectrum
01685                         KMO_TRY_EXIT_IF_NULL(
01686                             lambda_x = kmo_create_lambda_vec(gd.l.dim, 1,
01687                                                             gd.l.start,
01688                                                             gd.l.delta));
01689                         //
01690                         // spectrum correction
01691                         //
01692                         if ((strcmp(star_type, "O") == 0) ||
01693                             (strcmp(star_type, "B") == 0) ||
01694                             (strcmp(star_type, "A") == 0) ||
01695                             (strcmp(star_type, "F") == 0))
01696                         {
01697                             // we have a OBAF star
01698 
01699                             // if ATMOS_MODEL is present, lines will be removed
01700                             if (cpl_frameset_count_tags(frameset, ATMOS_MODEL) == 1) {
01701                                 // interpolate ATMOS_MODEL to same scale as data
01702                                 KMO_TRY_EXIT_IF_NULL(
01703                                     tmp_frame = kmo_dfs_get_frame(frameset, ATMOS_MODEL));
01704 
01705                                 KMO_TRY_EXIT_IF_NULL(
01706                                     atmos_model = kmo_interpolate_vector_wcs(tmp_frame, lambda_x));
01707 cpl_vector *tmp_spec_data_orig = NULL;
01708 int plot_it = 0;
01709 if (plot_it) {
01710     // store original spectrum
01711     KMO_TRY_EXIT_IF_NULL(
01712         tmp_spec_data_orig = cpl_vector_duplicate(tmp_spec_data));
01713 }
01714                                 // remove band-specific lines
01715                                 if (strcmp(filter_id, "H") == 0) {
01716                                     for (l = 0; l < nr_lines_h; l++) {
01717                                         KMO_TRY_EXIT_IF_ERROR(
01718                                             kmo_remove_line(tmp_spec_data, lambda_x, atmos_model, lines_center_h[l], lines_width_h[l]));
01719                                     }
01720                                 } else if (strcmp(filter_id, "HK") == 0) {
01721                                     for (l = 0; l < nr_lines_hk; l++) {
01722                                         KMO_TRY_EXIT_IF_ERROR(
01723                                             kmo_remove_line(tmp_spec_data, lambda_x, atmos_model, lines_center_hk[l], lines_width_hk[l]));
01724                                     }
01725                                 } else if (strcmp(filter_id, "K") == 0) {
01726                                     for (l = 0; l < nr_lines_k; l++) {
01727                                         KMO_TRY_EXIT_IF_ERROR(
01728                                             kmo_remove_line(tmp_spec_data, lambda_x, atmos_model, lines_center_k[l], lines_width_k[l]));
01729                                     }
01730                                 } else if (strcmp(filter_id, "IZ") == 0) {
01731                                     for (l = 0; l < nr_lines_iz; l++) {
01732                                         KMO_TRY_EXIT_IF_ERROR(
01733                                             kmo_remove_line(tmp_spec_data, lambda_x, atmos_model, lines_center_iz[l], lines_width_iz[l]));
01734                                     }
01735                                 } else if (strcmp(filter_id, "YJ") == 0) {
01736                                     for (l = 0; l < nr_lines_yj; l++) {
01737                                         KMO_TRY_EXIT_IF_ERROR(
01738                                             kmo_remove_line(tmp_spec_data, lambda_x, atmos_model, lines_center_yj[l], lines_width_yj[l]));
01739                                     }
01740                                 }
01741 if (plot_it) {
01742     cpl_vector *tmp_spec_data_atmo = NULL;
01743     cpl_vector *tmp_spec_data_new = NULL;
01744     KMO_TRY_EXIT_IF_NULL(
01745         tmp_spec_data_atmo = cpl_vector_duplicate(tmp_spec_data_orig));
01746     KMO_TRY_EXIT_IF_NULL(
01747         tmp_spec_data_new = cpl_vector_duplicate(tmp_spec_data));
01748     KMO_TRY_EXIT_IF_ERROR(
01749         cpl_vector_divide(tmp_spec_data_atmo, atmos_model));
01750 
01751     char *sss = cpl_sprintf("atmo_div_%s.fits", filter_id);
01752     if (i == 1) {
01753         cpl_vector_save(tmp_spec_data_atmo, sss, CPL_BPP_IEEE_DOUBLE, stored_sub_tel_data_headers[ifu_nr-1], CPL_IO_CREATE);
01754     } else {
01755         cpl_vector_save(tmp_spec_data_atmo, sss, CPL_BPP_IEEE_DOUBLE, stored_sub_tel_data_headers[ifu_nr-1], CPL_IO_EXTEND);
01756     }
01757 
01758     cpl_vector *med_vec = cpl_vector_duplicate(tmp_spec_data_orig);
01759     double  median = cpl_vector_get_median(med_vec);
01760     cpl_vector_delete(med_vec);
01761     int ii = 0;
01762     for (ii = 0; ii < cpl_vector_get_size(tmp_spec_data_orig); ii++) {
01763         if (cpl_vector_get(tmp_spec_data_orig, ii) < median/8)
01764             cpl_vector_set(tmp_spec_data_orig, ii, 0);
01765         if (cpl_vector_get(tmp_spec_data_atmo, ii) < median/8)
01766             cpl_vector_set(tmp_spec_data_atmo, ii, 0);
01767         if (cpl_vector_get(tmp_spec_data_new, ii) < median/8)
01768             cpl_vector_set(tmp_spec_data_new, ii, 0);
01769 
01770         if (cpl_vector_get(tmp_spec_data_orig, ii) > 3*median)
01771             cpl_vector_set(tmp_spec_data_orig, ii, 3*median);
01772         if (cpl_vector_get(tmp_spec_data_atmo, ii) > 3*median)
01773             cpl_vector_set(tmp_spec_data_atmo, ii, 3*median);
01774         if (cpl_vector_get(tmp_spec_data_new, ii) > 3*median)
01775             cpl_vector_set(tmp_spec_data_new, ii, 3*median);
01776     }
01777 
01778     double *pspec_dup = cpl_vector_get_data(tmp_spec_data_atmo);
01779     for (ii = 0; ii < cpl_vector_get_size(tmp_spec_data_atmo); ii++) {
01780         if (kmclipm_is_nan_or_inf(pspec_dup[ii])) {
01781             pspec_dup[ii] = 0.;
01782         }
01783     }
01784 
01785     cpl_bivector *plots[3];
01786     plots[0] = cpl_bivector_wrap_vectors((cpl_vector*)lambda_x, tmp_spec_data_orig);
01787     plots[1] = cpl_bivector_wrap_vectors((cpl_vector*)lambda_x, tmp_spec_data_atmo);
01788     plots[2] = cpl_bivector_wrap_vectors((cpl_vector*)lambda_x, tmp_spec_data_new);
01789     char *options[3] = {"w l t 'original'",
01790                         "w l t 'atmo divided'",
01791                         "w l t 'lines removed'"};
01792     sss = cpl_sprintf("set title '%s-band line removal (DET #%d)';", filter_id, i);
01793     cpl_plot_bivectors(sss,
01794                        (const char**)options, "", (const cpl_bivector**)plots, 3);
01795 //    cpl_plot_bivectors("set title 'Spectrum with lines removed'; set xrange [2.14:2.19];",
01796 //                       (const char**)options, "", (const cpl_bivector**)plots, 2);
01797     cpl_bivector_unwrap_vectors(plots[0]);
01798     cpl_bivector_unwrap_vectors(plots[1]);
01799     cpl_bivector_unwrap_vectors(plots[2]);
01800     cpl_free(sss); sss = NULL;
01801     cpl_vector_delete(tmp_spec_data_orig); tmp_spec_data_orig = NULL;
01802     cpl_vector_delete(tmp_spec_data_atmo); tmp_spec_data_atmo = NULL;
01803     cpl_vector_delete(tmp_spec_data_new); tmp_spec_data_new = NULL;
01804 }
01805                                 cpl_vector_delete(atmos_model); atmos_model = NULL;
01806                             } else {
01807                                 if (line_warning == FALSE) {
01808                                     cpl_msg_warning("", "No atmospheric model (ATMOS_MODEL) provided! "
01809                                                         "Won't remove any lines.");
01810                                     line_warning = TRUE;
01811                                 }
01812                             }
01813                         } else if (strcmp(star_type, "G") == 0) {
01814                             // we have a G star
01815                             if (cpl_frameset_count_tags(frameset, SOLAR_SPEC) == 1) {
01816                                 // interpolate SOLAR_SPEC to same scale as data
01817                                 // and divide it
01818                                 KMO_TRY_EXIT_IF_NULL(
01819                                     tmp_frame = kmo_dfs_get_frame(frameset, SOLAR_SPEC));
01820 
01821                                 // check if SOLAR_SPEC is the filter_id-one
01822                                 KMO_TRY_EXIT_IF_NULL(
01823                                     tmp_sub_header = kmclipm_propertylist_load(cpl_frame_get_filename(tmp_frame), 0));
01824                                 KMO_TRY_EXIT_IF_NULL(
01825                                     tmp_str = cpl_propertylist_get_string(tmp_sub_header, FILT_ID));
01826                                 KMO_TRY_ASSURE(strcmp(filter_id, tmp_str) == 0,
01827                                                CPL_ERROR_ILLEGAL_INPUT,
01828                                                "SOLAR_SPEC model must have primary "
01829                                                "keyword '%s' equal '%s'!!!",
01830                                                FILT_ID, filter_id);
01831                                 cpl_propertylist_delete(tmp_sub_header); tmp_sub_header = NULL;
01832 
01833                                 KMO_TRY_EXIT_IF_NULL(
01834                                     solar_spec = kmo_interpolate_vector_wcs(tmp_frame, lambda_x));
01835 
01836                                 // values are set to zero if solar_spec isn't
01837                                 // overlapping wavelength range of star apectrum
01838                                 // completely
01839                                 KMO_TRY_EXIT_IF_ERROR(
01840                                     cpl_vector_divide(tmp_spec_data, solar_spec));
01841                                 cpl_vector_delete(solar_spec); solar_spec = NULL;
01842                             } else {
01843                                 if (print_warning_once == TRUE) {
01844                                     cpl_msg_warning("","No solar spectrum (SOLAR_SPEC) provided! "
01845                                                        "Can't divide it from extracted "
01846                                                        "standard star spectrum!");
01847                                     print_warning_once = FALSE;
01848                                 }
01849                             }
01850                         } else {
01851 //                            cpl_msg_warning("","No startype was provided! Can't"
01852 //                                            " divide solar spectrum for G stars "
01853 //                                            "or fit a profile to atmospheric "
01854 //                                            "transmission for OBAF stars.");
01855                         }
01856 
01857                         if (star_temperature > 0.0) {
01858                             // divide blackbody from tmp_spec_data
01859                             KMO_TRY_EXIT_IF_ERROR(
01860                                 kmo_divide_blackbody(tmp_spec_data, lambda_x, star_temperature));
01861                         }
01862 
01863                         cpl_vector_delete(lambda_x); lambda_x = NULL;
01864 
01865                         // normalise telluric and its noise
01866                         // mean is taken in lambda defined range
01867                         KMO_TRY_EXIT_IF_ERROR(
01868                             kmo_calc_band_mean(stored_sub_tel_data_headers[ifu_nr-1],
01869                                                filter_id,
01870                                                tmp_spec_data,
01871                                                tmp_spec_noise,
01872                                                &mean_data,
01873                                                &mean_noise));
01874 
01875                         KMO_TRY_EXIT_IF_ERROR(
01876                             cpl_vector_divide_scalar(tmp_spec_data, mean_data));
01877 
01878                         if (tmp_spec_noise != NULL) {
01879                             KMO_TRY_EXIT_IF_ERROR(
01880                                 cpl_vector_divide_scalar(tmp_spec_noise, mean_noise));
01881 
01882                             // set noise spectrum also to zero when solar_spec is too short
01883                             KMO_TRY_EXIT_IF_NULL(
01884                                 ptmp_spec_data = cpl_vector_get_data_const(tmp_spec_data));
01885                             KMO_TRY_EXIT_IF_NULL(
01886                                 ptmp_spec_noise = cpl_vector_get_data(tmp_spec_noise));
01887                             for (k = 0; k < cpl_vector_get_size(tmp_spec_data); k++) {
01888                                 if (ptmp_spec_data[k] == 0.0) {
01889                                     ptmp_spec_noise[k] = 0.0;
01890                                 }
01891                             }
01892                         }
01893                         KMO_TRY_CHECK_ERROR_STATE();
01894 
01895                         // store telluric & error spectrum
01896                         stored_telluric_data[ifu_nr-1] = tmp_spec_data;
01897                         stored_telluric_noise[ifu_nr-1] = tmp_spec_noise;
01898 
01899                         // if magnitude is provided
01900                         // calculate zeropoint and throughput
01901                         if (has_magnitude) {
01902                             // calculate QC THROUGHPUT
01903                             crpix1 = cpl_propertylist_get_double(stored_sub_tel_data_headers[ifu_nr-1], CRPIX1);
01904                             crval1 = cpl_propertylist_get_double(stored_sub_tel_data_headers[ifu_nr-1], CRVAL1);
01905                             cdelt1 = cpl_propertylist_get_double(stored_sub_tel_data_headers[ifu_nr-1], CDELT1);
01906                             KMO_TRY_CHECK_ERROR_STATE();
01907 
01908                             KMO_TRY_EXIT_IF_ERROR(
01909                                 kmo_calc_counts(spec_qc, filter_id,
01910                                                 crpix1, crval1, cdelt1,
01911                                                 &counts1, &counts2));
01912                             KMO_TRY_CHECK_ERROR_STATE();
01913 
01914                             counts1 /= exptime;
01915                             counts2 /= exptime;
01916 
01917                             stored_qc_throughput[ifu_nr-1] =
01918                                 kmo_calc_throughput(magnitude1, magnitude2, counts1, counts2,
01919                                                     cpl_propertylist_get_double(stored_sub_tel_data_headers[ifu_nr-1], GAIN),
01920                                                     filter_id);
01921                             KMO_TRY_CHECK_ERROR_STATE();
01922 
01923                             if (kmclipm_is_nan_or_inf(stored_qc_throughput[ifu_nr-1])) {
01924                                 stored_qc_throughput[ifu_nr-1] = -1;
01925                             }
01926                             KMO_TRY_EXIT_IF_ERROR(
01927                                 kmclipm_update_property_double(stored_sub_tel_data_headers[ifu_nr-1],
01928                                                                QC_THROUGHPUT,
01929                                                                stored_qc_throughput[ifu_nr-1],
01930                                                                "[] IFU throughput"));
01931 
01932                             // calculate QC ZEROPOINT
01933                             zeropoint = kmo_calc_zeropoint(magnitude1, magnitude2, counts1, counts2, cdelt3, filter_id);
01934                             if (kmclipm_is_nan_or_inf(zeropoint)) {
01935                                 zeropoint = -1;
01936                             }
01937                             KMO_TRY_CHECK_ERROR_STATE();
01938 
01939                             KMO_TRY_EXIT_IF_ERROR(
01940                                 kmclipm_update_property_double(stored_sub_tel_data_headers[ifu_nr-1],
01941                                                                QC_ZEROPOINT,
01942                                                                zeropoint,
01943                                                                "[mag] IFU zeropoint"));
01944                         }
01945                         cpl_vector_delete(spec_qc);spec_qc = NULL;
01946                     } else {
01947                         cpl_error_reset();
01948                         // IFU is invalid
01949                         KMO_TRY_EXIT_IF_NULL(
01950                             stored_sub_tel_data_headers[ifu_nr-1] =
01951                                 cpl_propertylist_duplicate(sub_header_orig));
01952                         KMO_TRY_EXIT_IF_NULL(
01953                             stored_sub_tel_noise_headers[ifu_nr-1] =
01954                                 cpl_propertylist_duplicate(sub_header_orig));
01955                         KMO_TRY_EXIT_IF_NULL(
01956                             stored_sub_psf_headers[ifu_nr-1] =
01957                                 cpl_propertylist_duplicate(sub_header_orig));
01958                         if (save_cubes) {
01959                             KMO_TRY_EXIT_IF_NULL(
01960                                 stored_sub_cube_data_headers[ifu_nr-1] =
01961                                     cpl_propertylist_duplicate(sub_header_orig));
01962                             KMO_TRY_EXIT_IF_NULL(
01963                                 stored_sub_cube_noise_headers[ifu_nr-1] =
01964                                     cpl_propertylist_duplicate(sub_header_orig));
01965                         }
01966                     }
01967 
01968                     // create EXTNAME keyword as DATA
01969                     KMO_TRY_EXIT_IF_NULL(
01970                         extname = kmo_extname_creator(ifu_frame, ifu_nr, EXT_DATA));
01971                     KMO_TRY_EXIT_IF_ERROR(
01972                         kmclipm_update_property_string(stored_sub_tel_data_headers[ifu_nr-1],
01973                                                        EXTNAME, extname, "FITS extension name"));
01974                     KMO_TRY_EXIT_IF_ERROR(
01975                         kmclipm_update_property_string(stored_sub_psf_headers[ifu_nr-1],
01976                                                        EXTNAME, extname, "FITS extension name"));
01977                     if (save_cubes) {
01978                         KMO_TRY_EXIT_IF_ERROR(
01979                             kmclipm_update_property_string(stored_sub_cube_data_headers[ifu_nr-1],
01980                                                            EXTNAME, extname, "FITS extension name"));
01981                     }
01982                     cpl_free(extname); extname = NULL;
01983 
01984                     // create EXTNAME keyword as NOISE
01985                     if (stored_sub_tel_noise_headers[ifu_nr-1] == NULL) {
01986                         KMO_TRY_EXIT_IF_NULL(
01987                             stored_sub_tel_noise_headers[ifu_nr-1] =
01988                                     cpl_propertylist_duplicate(
01989                                             stored_sub_tel_data_headers[ifu_nr-1]));
01990                     }
01991                     KMO_TRY_EXIT_IF_NULL(
01992                         extname = kmo_extname_creator(ifu_frame, ifu_nr, EXT_NOISE));
01993                     KMO_TRY_EXIT_IF_ERROR(
01994                         kmclipm_update_property_string(stored_sub_tel_noise_headers[ifu_nr-1],
01995                                                        EXTNAME, extname, "FITS extension name"));
01996                     if (save_cubes) {
01997                         KMO_TRY_EXIT_IF_NULL(
01998                             stored_sub_cube_noise_headers[ifu_nr-1] =
01999                                     cpl_propertylist_duplicate(
02000                                             stored_sub_cube_data_headers[ifu_nr-1]));
02001                         KMO_TRY_EXIT_IF_ERROR(
02002                             kmclipm_update_property_string(stored_sub_cube_noise_headers[ifu_nr-1],
02003                                                            EXTNAME, extname, "FITS extension name"));
02004                     }
02005                     cpl_free(extname); extname = NULL;
02006                 } // for j ifus (load, process & store)
02007                 cpl_propertylist_delete(sub_header_orig); sub_header_orig = NULL;
02008             } // for i detectors (load, process & store)
02009             KMO_TRY_CHECK_ERROR_STATE();
02010 
02011             // write QC parameter: nr of std stars
02012             KMO_TRY_EXIT_IF_ERROR(
02013                 kmclipm_update_property_int(main_header_tel, QC_NR_STD_STARS,
02014                                         nr_std_stars, "[] Nr. of std stars"));
02015 
02016             // update which IFUs are not used
02017             kmo_print_unused_ifus(unused_ifus_after, TRUE);
02018 
02019             KMO_TRY_EXIT_IF_ERROR(
02020                 kmo_set_unused_ifus(unused_ifus_after, main_header_tel, "kmo_std_star"));
02021 
02022             KMO_TRY_EXIT_IF_NULL(
02023                 main_header_psf = cpl_propertylist_duplicate(main_header_tel));
02024 
02025             if (has_magnitude) {
02026                 // calculate QC THROUGHPUT MEAN and QC THROUGHPUT SDV
02027                 // and update main header
02028                 KMO_TRY_EXIT_IF_ERROR(
02029                     kmo_calc_mean_throughput(stored_qc_throughput,
02030                                              nr_devices * KMOS_IFUS_PER_DETECTOR,
02031                                              &throughput_mean,
02032                                              &throughput_sdv));
02033                 KMO_TRY_EXIT_IF_ERROR(
02034                     kmclipm_update_property_double(main_header_tel, QC_THROUGHPUT_MEAN,
02035                                                    throughput_mean, "[] mean throughput for all detectors"));
02036                 KMO_TRY_EXIT_IF_ERROR(
02037                     kmclipm_update_property_double(main_header_tel, QC_THROUGHPUT_SDV,
02038                                                    throughput_sdv, "[] stdev throughput for all detectors"));
02039             }
02040             KMO_TRY_CHECK_ERROR_STATE();
02041 
02042             //
02043             // save output data
02044             //
02045             if (!suppress_extension) {
02046                 KMO_TRY_EXIT_IF_NULL(
02047                     fn_suffix = cpl_sprintf("%s", suffix));
02048             } else {
02049                 KMO_TRY_EXIT_IF_NULL(
02050                     fn_suffix = cpl_sprintf("%s", ""));
02051             }
02052 
02053             // save primary extension
02054             cpl_msg_info("","Saving STD exposure No. %d", nr_exp+1);
02055             KMO_TRY_EXIT_IF_ERROR(
02056                 kmo_dfs_save_main_header(frameset, filename_telluric, fn_suffix,
02057                                          obj_frame, main_header_tel, parlist,
02058                                          cpl_func));
02059             KMO_TRY_EXIT_IF_ERROR(
02060                 kmo_dfs_save_main_header(frameset, filename_starspec, fn_suffix,
02061                                          obj_frame, main_header_tel, parlist,
02062                                          cpl_func));
02063             KMO_TRY_EXIT_IF_ERROR(
02064                 kmo_dfs_save_main_header(frameset, filename_mask, fn_suffix,
02065                                          obj_frame, main_header_psf, parlist,
02066                                          cpl_func));
02067             KMO_TRY_EXIT_IF_ERROR(
02068                 kmo_dfs_save_main_header(frameset, filename_psf, fn_suffix,
02069                                          obj_frame, main_header_psf, parlist,
02070                                          cpl_func));
02071             if (save_cubes) {
02072                 KMO_TRY_EXIT_IF_ERROR(
02073                     kmo_dfs_save_main_header(frameset, filename_cubes, fn_suffix,
02074                                              obj_frame, main_header_psf, parlist,
02075                                              cpl_func));
02076             }
02077 
02078             // save stored frames
02079             for (i = 1; i <= nr_devices; i++) {
02080                 for (j = 0; j < KMOS_IFUS_PER_DETECTOR; j++) {
02081                     ifu_nr = (i-1)*KMOS_IFUS_PER_DETECTOR + j + 1;
02082 
02083                     // save telluric-vector
02084 kmclipm_vector *ddd = NULL;
02085 if (stored_telluric_data[ifu_nr-1] != NULL)
02086 ddd = kmclipm_vector_create(cpl_vector_duplicate(stored_telluric_data[ifu_nr-1]));
02087                     KMO_TRY_EXIT_IF_ERROR(
02088                         kmo_dfs_save_vector(ddd, filename_telluric, fn_suffix,
02089                                             stored_sub_tel_data_headers[ifu_nr-1],
02090                                             0./0.));
02091 kmclipm_vector_delete(ddd); ddd =NULL;
02092 
02093 if (stored_telluric_noise[ifu_nr-1] != NULL)
02094 ddd = kmclipm_vector_create(cpl_vector_duplicate(stored_telluric_noise[ifu_nr-1]));
02095                     KMO_TRY_EXIT_IF_ERROR(
02096                         kmo_dfs_save_vector(ddd, filename_telluric, fn_suffix,
02097                                             stored_sub_tel_noise_headers[ifu_nr-1],
02098                                             0./0.));
02099 kmclipm_vector_delete(ddd); ddd =NULL;
02100 
02101                     // save star_spec-vector
02102 if (stored_starspec_data[ifu_nr-1] != NULL)
02103 ddd = kmclipm_vector_create(cpl_vector_duplicate(stored_starspec_data[ifu_nr-1]));
02104                     KMO_TRY_EXIT_IF_ERROR(
02105                         kmo_dfs_save_vector(ddd, filename_starspec, fn_suffix,
02106                                             stored_sub_tel_data_headers[ifu_nr-1],
02107                                             0./0.));
02108 kmclipm_vector_delete(ddd); ddd =NULL;
02109 
02110 if (stored_starspec_noise[ifu_nr-1] != NULL)
02111 ddd = kmclipm_vector_create(cpl_vector_duplicate(stored_starspec_noise[ifu_nr-1]));
02112                     KMO_TRY_EXIT_IF_ERROR(
02113                         kmo_dfs_save_vector(ddd, filename_starspec, fn_suffix,
02114                                             stored_sub_tel_noise_headers[ifu_nr-1],
02115                                             0./0.));
02116 kmclipm_vector_delete(ddd); ddd =NULL;
02117 
02118                     // save psf-image
02119                     KMO_TRY_EXIT_IF_ERROR(
02120                         kmo_dfs_save_image(stored_psf_data[ifu_nr-1],
02121                                            filename_psf, fn_suffix,
02122                                            stored_sub_psf_headers[ifu_nr-1],
02123                                            0./0.));
02124 
02125                     // save mask-image
02126                     KMO_TRY_EXIT_IF_ERROR(
02127                         kmo_dfs_save_image(stored_mask[ifu_nr-1],
02128                                            filename_mask, fn_suffix,
02129                                            stored_sub_psf_headers[ifu_nr-1],
02130                                            0./0.));
02131                     // save reonstructed cubes
02132                     if (save_cubes) {
02133                         KMO_TRY_EXIT_IF_ERROR(
02134                             kmo_dfs_save_cube(stored_data_cube[ifu_nr-1],
02135                                               filename_cubes, fn_suffix,
02136                                               stored_sub_cube_data_headers[ifu_nr-1],
02137                                               0./0.));
02138                         KMO_TRY_EXIT_IF_ERROR(
02139                             kmo_dfs_save_cube(stored_noise_cube[ifu_nr-1],
02140                                               filename_cubes, fn_suffix,
02141                                               stored_sub_cube_noise_headers[ifu_nr-1],
02142                                               0./0.));
02143                     }
02144                 } // for j ifus (save stored)
02145             } // for i detectors (save stored)
02146             KMO_TRY_CHECK_ERROR_STATE();
02147         } // if (frameCnt == 0)
02148     }
02149     KMO_CATCH
02150     {
02151         KMO_CATCH_MSG();
02152         ret_val = -1;
02153     }
02154 
02155     kmo_delete_objSkyStruct(obj_sky_struct);
02156     kmo_free_fits_desc(&desc1);
02157     kmo_free_fits_desc(&desc2);
02158     kmo_free_unused_ifus(unused_ifus_before); unused_ifus_before = NULL;
02159     kmo_free_unused_ifus(unused_ifus_after); unused_ifus_after = NULL;
02160     cpl_free(bounds); bounds = NULL;
02161     cpl_propertylist_delete(main_header_tel); main_header_tel = NULL;
02162     cpl_propertylist_delete(main_header_psf); main_header_psf = NULL;
02163     cpl_vector_delete(atmos_model); atmos_model = NULL;
02164     cpl_vector_delete(solar_spec); solar_spec = NULL;
02165     cpl_table_delete(spec_type_LUT); spec_type_LUT = NULL;
02166     cpl_vector_delete(identified_slices); identified_slices = NULL;
02167     cpl_propertylist_delete(sub_header_orig); sub_header_orig = NULL;
02168     for (i = 0; i < nr_devices * KMOS_IFUS_PER_DETECTOR; i++) {
02169         cpl_vector_delete(stored_telluric_data[i]); stored_telluric_data[i] = NULL;
02170         cpl_vector_delete(stored_telluric_noise[i]); stored_telluric_noise[i] = NULL;
02171         cpl_vector_delete(stored_starspec_data[i]); stored_starspec_data[i] = NULL;
02172         cpl_vector_delete(stored_starspec_noise[i]); stored_starspec_noise[i] = NULL;
02173         cpl_image_delete(stored_psf_data[i]); stored_psf_data[i] = NULL;
02174         cpl_propertylist_delete(stored_sub_tel_data_headers[i]); stored_sub_tel_data_headers[i] = NULL;
02175         cpl_propertylist_delete(stored_sub_tel_noise_headers[i]); stored_sub_tel_noise_headers[i] = NULL;
02176         if (save_cubes) {
02177             cpl_propertylist_delete(stored_sub_cube_data_headers[i]); stored_sub_cube_data_headers[i] = NULL;
02178             cpl_propertylist_delete(stored_sub_cube_noise_headers[i]); stored_sub_cube_noise_headers[i] = NULL;
02179         }
02180         cpl_propertylist_delete(stored_sub_psf_headers[i]); stored_sub_psf_headers[i] = NULL;
02181         cpl_image_delete(stored_mask[i]); stored_mask[i] = NULL;
02182         cpl_imagelist_delete(stored_data_cube[i]); stored_data_cube[i] = NULL;
02183         cpl_imagelist_delete(stored_noise_cube[i]); stored_noise_cube[i] = NULL;
02184     }
02185     cpl_free(stored_telluric_data); stored_telluric_data = NULL;
02186     cpl_free(stored_telluric_noise); stored_telluric_noise = NULL;
02187     cpl_free(stored_starspec_data); stored_starspec_data = NULL;
02188     cpl_free(stored_starspec_noise); stored_starspec_noise = NULL;
02189     cpl_free(stored_psf_data); stored_psf_data = NULL;
02190     cpl_free(stored_sub_tel_data_headers); stored_sub_tel_data_headers = NULL;
02191     cpl_free(stored_sub_tel_noise_headers); stored_sub_tel_noise_headers = NULL;
02192     if (save_cubes) {
02193         cpl_free(stored_sub_cube_data_headers); stored_sub_cube_data_headers = NULL;
02194         cpl_free(stored_sub_cube_noise_headers); stored_sub_cube_noise_headers = NULL;
02195     }
02196     cpl_free(stored_sub_psf_headers); stored_sub_psf_headers = NULL;
02197     cpl_free(stored_qc_throughput); stored_qc_throughput = NULL;
02198     cpl_free(suffix); suffix = NULL;
02199     cpl_free(fn_suffix); fn_suffix = NULL;
02200     cpl_free(stored_mask); stored_mask = NULL;
02201     cpl_free(stored_data_cube); stored_data_cube = NULL;
02202     cpl_free(stored_noise_cube); stored_noise_cube = NULL;
02203     cpl_free(grat_id); grat_id = NULL;
02204     cpl_frameset_delete(frameset_std); frameset_std = NULL;
02205 
02206     return ret_val;
02207 }
02208