sinfo_new_wave_cal_slit2.c

00001 /*
00002  * This file is part of the ESO SINFONI Pipeline
00003  * Copyright (C) 2004,2005 European Southern Observatory
00004  *
00005  * This program is free software; you can redistribute it and/or modify
00006  * it under the terms of the GNU General Public License as published by
00007  * the Free Software Foundation; either version 2 of the License, or
00008  * (at your option) any later version.
00009  *
00010  * This program is distributed in the hope that it will be useful,
00011  * but WITHOUT ANY WARRANTY; without even the implied warranty of
00012  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00013  * GNU General Public License for more details.
00014  *
00015  * You should have received a copy of the GNU General Public License
00016  * along with this program; if not, write to the Free Software
00017  * Foundation, 51 Franklin St, Fifth Floor, Boston, MA  02111-1307  USA
00018  */
00019 /*----------------------------------------------------------------------------
00020  
00021    File name    :       sinfo_new_wave_cal_slit2.c
00022    Author       :    A. Modigliani
00023    Created on   :    Sep 17, 2003
00024    Description  : 
00025 
00026  ---------------------------------------------------------------------------*/
00027 
00028 #ifdef HAVE_CONFIG_H
00029 #  include <config.h>
00030 #endif
00031 /*----------------------------------------------------------------------------
00032                                 Includes
00033  ---------------------------------------------------------------------------*/
00034 #include "sinfo_new_wave_cal_slit2.h"
00035 #include "sinfo_pro_save.h"
00036 #include "sinfo_pro_types.h"
00037 #include "sinfo_wavecal_ini_by_cpl.h"
00038 #include "sinfo_wcal_functions.h"
00039 #include "sinfo_absolute.h"
00040 #include "sinfo_wave_calibration.h"
00041 #include "sinfo_wavecal.h"
00042 #include "sinfo_globals.h"
00043 #include "sinfo_hidden.h"
00044 
00045 #include "sinfo_utilities.h"
00046 #include "sinfo_utils_wrappers.h"
00047 #include "sinfo_error.h"
00048 
00049 
00050 /*----------------------------------------------------------------------------
00051                                 Defines
00052  ---------------------------------------------------------------------------*/
00053 
00054 /*----------------------------------------------------------------------------
00055                              Function Definitions
00056  ---------------------------------------------------------------------------*/
00057 static cpl_error_code
00058 sinfo_image_resample(const char* plugin_id,
00059                      cpl_parameterlist* config,
00060                      cpl_frameset* sof,
00061                      cpl_frameset* ref_set);
00062 
00070 /*----------------------------------------------------------------------------
00071    Function     :       sinfo_wave_cal_slit()
00072    In           :       ini_file: file name of according .ini file
00073    Out          :       integer (0 if it worked, -1 if it doesn't) 
00074    Job          :
00075 
00076 
00077   Normal method:
00078 
00079   does the wavelength calibration and the fitting of the slitlet sinfo_edge 
00080   positions (ASCII file 32 x 2 values) if wished
00081   produces an array of the bcoefs and of the fit parameters if wished and a 
00082   wavelength calibration map input is an emission line frame and a line list
00083 
00084 
00085   o searching for lines by cross sinfo_correlation with a line list
00086   o Gaussian fitting of emission lines in each column->positions of the lines->
00087     resulting fit parameters can be stored in an ASCII file
00088   o Fitting of a polynomial to the line positions for each column
00089   o Smoothing: fitting of each polynomial coefficient by another polynomial
00090     across the whole frame -> resulting polynomial coefficients can be stored 
00091     in an ASCII file.
00092   o Wavelength calibration map (micron value for each frame pixel) can be
00093     produced by using these coefficients and a cross sinfo_correlation to the
00094     original frame
00095 
00096   o The slitlet sinfo_edge positions can be fitted:
00097     1) Automatically (not really stable) or by using guess sinfo_edge positions
00098     2) By using a Boltzmann or a linear slope function
00099 
00100   Slit method:
00101 
00102   does the wavelength calibration and the fitting of the slitlet sinfo_edge 
00103   positions (ASCII file 32 x 2 values) if wished produces a list of the fit 
00104   parameters and of the smoothed coefficients if wished and a wavelength 
00105   calibration map input is an emission line frame and a line list
00106 
00107   o Does the same as other method but smoothes the found polynomial 
00108     coefficients within each slitlet and not over the whole frame.
00109 
00110   o Produces always a wavelength calibration map and does not crosscorrelate.
00111 
00112  ---------------------------------------------------------------------------*/
00113 
00114 
00115 int 
00116 sinfo_new_wave_cal_slit2(const char* plugin_id,
00117                         cpl_parameterlist* config, 
00118                         cpl_frameset* sof,cpl_frameset* ref_set)
00119 {
00120   wave_config * cfg =NULL;
00121   char col_name[MAX_NAME_SIZE];
00122   char tbl_name[MAX_NAME_SIZE];
00123   char tbl_fitpar_name[MAX_NAME_SIZE];
00124   char tbl_line_list_name[MAX_NAME_SIZE];
00125   char tbl_slitpos_guess_name[MAX_NAME_SIZE];
00126   char key_name[MAX_NAME_SIZE];
00127   char col[MAX_NAME_SIZE];
00128   int check = 0;
00129   int lx = 0;
00130   int ly = 0;
00131   int n_lines=0;
00132   int i = 0;
00133   int j = 0;
00134   int n = 0;
00135   int readsum=0;
00136   int sum=0;
00137   int fit=0;
00138   int sw=0;
00139  
00140   int* status=NULL;
00141   int* n_found_lines=NULL;
00142   int* sum_pointer=NULL;
00143   int** row_clean=NULL;
00144 
00145   float a=0;
00146   float shift=0;
00147   float val_x=0;
00148   float val_y=0;
00149 
00150   float* wave=NULL;
00151   float* intens=NULL;
00152 
00153   float** acoefs=NULL;
00154   float** wavelength_clean=NULL;
00155   float** sinfo_slit_pos=NULL;
00156 
00157   double fwhm_med=0;
00158   double fwhm_avg=0;
00159   double coef_med=0;
00160   double coef_avg=0;
00161  
00162   cpl_image * im=NULL ;
00163 
00164   FitParams** par=NULL;
00165 
00166   cpl_table* tbl_wcal=NULL;
00167   cpl_table* tbl_spos=NULL;
00168   cpl_table* tbl_fitpar = NULL;
00169   cpl_table* tbl_line_list = NULL;
00170   cpl_table* tbl_slitpos_guess=NULL;
00171   cpl_table * tbl_fp =NULL;
00172   cpl_table* qclog_tbl=NULL;
00173 
00174   cpl_image* map_img=NULL;
00175 
00176   cpl_frameset* raw=NULL;
00177   cpl_parameter* p=NULL;
00178 
00179   qc_wcal* qc=sinfo_qc_wcal_new();
00180   int pdensity=0;
00181 
00182   /*   -----------------------------------------------------------------
00183        1) parse the file names and parameters to the ns_config data 
00184           structure cfg
00185        -----------------------------------------------------------------
00186   */
00187 
00188   check_nomsg(p=cpl_parameterlist_find(config,"sinfoni.product.density"));
00189   check_nomsg(pdensity=cpl_parameter_get_int(p));
00190 
00191 
00192   sinfo_msg("Parsing cpl input");
00193   check_nomsg(raw=cpl_frameset_new());
00194   cknull(cfg = sinfo_parse_cpl_input_wave(config,sof,&raw),
00195      "could not parse cpl input!") ;
00196  
00197   check_nomsg(p = cpl_parameterlist_find(config,
00198                                          "sinfoni.wavecal.slitpos_boostrap"));
00199   check_nomsg(sw=cpl_parameter_get_bool(p));
00200 
00201   if(sw==1) {
00202      cfg->nslitlets=32;
00203      cfg->calibIndicator=1;
00204      cfg->wavemapInd=0;
00205      cfg->slitposIndicator=1;
00206      sinfo_msg("***********************************");
00207      sinfo_msg("parameter setting for %s",PRO_WAVE_SLITPOS_STACKED);
00208      sinfo_msg("***********************************");
00209   }
00210 
00211   if(sinfo_is_fits_file(cfg->inFrame) != 1) {
00212     sinfo_msg_error("Input file cfg->inFrame %s is not FITS",cfg->inFrame);
00213     goto cleanup;
00214   }
00215 
00216 
00217   if (cfg->slitposIndicator == 1 && cfg->estimateIndicator == 1) {
00218     if (sinfo_is_fits_file(cfg->slitposGuessName) != 1) {
00219       sinfo_msg_error("slitlet position guess list not given!");
00220       goto cleanup;
00221     }
00222   }
00223 
00224   if (cfg->calibIndicator == 0 && cfg->wavemapInd == 1) {
00225     if (sinfo_is_fits_file(cfg->coeffsName) != 1) {
00226       sinfo_msg_error("coefficients list not given!");
00227       goto cleanup;
00228     }
00229   }
00230 
00231   if (cfg->slitposIndicator == 1) {
00232     if (cfg->calibIndicator != 1 && cfg->estimateIndicator != 1) {
00233       if (sinfo_is_fits_file(cfg->paramsList) != 1) {
00234     sinfo_msg_error("parameter list not given!");
00235     goto cleanup;
00236       }
00237     }
00238   }
00239  
00240   /*---load the emission line frame---*/
00241   check(im = cpl_image_load(cfg->inFrame,CPL_TYPE_FLOAT,0,0)
00242         ,"could not load image");
00243   lx = cpl_image_get_size_x(im);
00244   ly = cpl_image_get_size_y(im);
00245 
00246   if (cfg->calibIndicator == 1 || cfg->wavemapInd == 1) {
00247     /*---open the line list and read the number of lines---*/
00248     strcpy(tbl_line_list_name,cfg->lineList);
00249     check_nomsg(tbl_line_list = cpl_table_load(tbl_line_list_name,1,0));
00250     check_nomsg(n = cpl_table_get_nrow(tbl_line_list));
00251     n_lines = n;
00252 
00253     check_nomsg(wave   = cpl_table_get_data_float(tbl_line_list,"wave"));
00254     check_nomsg(intens = cpl_table_get_data_float(tbl_line_list,"int"));
00255   }
00256 
00257   /*
00258   ----------------------------------------------------------------------
00259   ---------------------------FINDLINES----------------------------------
00260   ----------------------------------------------------------------------
00261   */
00262 
00263   /*if not yet done:
00264     do the wavelength calibration, that means: 
00265     find the dispersion relation and parameterize its coefficients
00266   */
00267   /*
00268   sinfo_msg("guessBeginWave=%g",cfg->guessBeginWavelength);
00269   sinfo_msg("guessDisp1=%g",cfg->guessDispersion1);
00270   sinfo_msg("guessDisp2=%g",cfg->guessDispersion2);
00271   */
00272   if (cfg->calibIndicator == 1 && cfg->wavemapInd == 0) {
00273     sinfo_msg("Findlines");
00274     acoefs  = sinfo_new_2Dfloatarray(cfg->nrDispCoefficients, lx);
00275 
00276     /*allocate memory*/
00277     n_found_lines    = sinfo_new_intarray(lx); 
00278     row_clean        = sinfo_new_2Dintarray(lx, n_lines);
00279     wavelength_clean = sinfo_new_2Dfloatarray(lx, n_lines);
00280     sum_pointer      = sinfo_new_intarray(1) ;
00281 
00282     /*find the emission lines in each image column*/
00283     sinfo_new_intarray_set_value(sum_pointer, 0, 0);
00284  
00285     ck0(check = sinfo_new_find_lines(im, 
00286                                 wave, 
00287                                 intens, 
00288                                 n_lines, 
00289                                 row_clean, 
00290                                 wavelength_clean, 
00291                                 cfg->guessBeginWavelength, 
00292                         cfg->guessDispersion1, 
00293                                 cfg->guessDispersion2,
00294                                 cfg->mindiff, 
00295                                 cfg->halfWidth, 
00296                                 n_found_lines, 
00297                                 cfg->sigma, 
00298                                 sum_pointer),
00299     "sinfo_findLines failed!");
00300 
00301     /*---------------------------------------------------------------------
00302      *-----------------------WAVE_CALIB-------------------------------------
00303      *---------------------------------------------------------------------
00304     */
00305     sinfo_msg("Wave Calibration");
00306     sum = sinfo_new_intarray_get_value(sum_pointer,0);
00307     /* allocate memory for the fit parameters */
00308     cknull(par = sinfo_new_fit_params( sum ),
00309        "sinfo_newFitParams failed!");
00310 
00311   /*
00312    fit each line, make a polynomial fit and fit the resulting fit 
00313    coefficients across the columns of the slitlet
00314    */
00315    cknull(map_img = sinfo_new_wave_cal(im, 
00316                   par, 
00317                   acoefs,
00318                   cfg->nslitlets, 
00319                   row_clean,
00320                   wavelength_clean,
00321                   n_found_lines,
00322                   cfg->guessDispersion1,
00323                   cfg->halfWidth,
00324                   cfg->minAmplitude,
00325                   cfg->maxResidual,
00326                   cfg->fwhm,
00327                   cfg->nrDispCoefficients, 
00328                   cfg->nrCoefCoefficients,
00329                   cfg->sigmaFactor,
00330                   cfg->pixeldist,
00331                   cfg->pixel_tolerance),
00332       "sinfo_wave_cal failed!");
00333 
00334    sinfo_msg("Check line positions");
00335 
00336    shift=sinfo_new_check_line_positions(im,acoefs,cfg->nrDispCoefficients,
00337                     cfg->guessDispersion1,par);
00338    if (FLAG == shift){ 
00339       sinfo_msg_error("checkForLinePositions failed!\n");
00340    }
00341 
00342 
00343    sinfo_det_ncounts(raw, cfg->qc_thresh_max, qc);
00344 
00345    cknull_nomsg(qclog_tbl = sinfo_qclog_init());
00346    ck0_nomsg(sinfo_qclog_add_int(qclog_tbl,"QC WAVE ALL",n_lines,
00347                                  "Number of found lines","%d"));
00348    ck0_nomsg(sinfo_qclog_add_int(qclog_tbl,"QC WAVE NPIXSAT",qc->nsat,
00349                                  "Number of saturated pixels","%d"));
00350    ck0_nomsg(sinfo_qclog_add_double(qclog_tbl,"QC WAVE MAXFLUX",qc->max_di,
00351                                  "Max int off-lamp subtracted frm","%g"));
00352    ck0_nomsg(sinfo_qclog_add_double(qclog_tbl,"QC WAVE POSERR",shift,
00353                                  "Overall positioning error in mum","%g"));
00354 
00355    ck0(sinfo_pro_save_ima(map_img,ref_set,sof,cfg->outName,
00356               PRO_WAVE_MAP,qclog_tbl,plugin_id,config),
00357        "cannot save ima %s", cfg->outName);
00358 
00359    sinfo_free_table(&qclog_tbl);
00360    sinfo_free_image(&map_img);
00361 
00362    /*
00363     #store the resulting polynomial fit coefficients in an 
00364      ASCII file if wished
00365    */
00366 
00367    if (cfg->writeCoeffsInd == 1) {
00368      check_nomsg(tbl_wcal = cpl_table_new(lx));
00369      for (i=0; i< cfg->nrDispCoefficients; i++) {
00370        snprintf(col_name,MAX_NAME_SIZE-1,"%s%d","coeff",i);
00371        check_nomsg(cpl_table_new_column(tbl_wcal,col_name, CPL_TYPE_DOUBLE));
00372      }
00373 
00374      cknull_nomsg(qclog_tbl = sinfo_qclog_init());
00375      ck0_nomsg(sinfo_qclog_add_int(qclog_tbl,"QC WAVE ALL",n_lines,
00376                                    "Number found lines","%d"));
00377 
00378      for (j=0; j< lx; j++) {
00379        for (i=0; i< cfg->nrDispCoefficients; i++) {
00380      snprintf(col_name,MAX_NAME_SIZE-1,"%s%d","coeff",i);
00381      a = sinfo_new_array2D_get_value(acoefs, i, j);
00382      /* fprintf(acoefs_file, "%15.13g ", a) ; */
00383      cpl_table_set_double(tbl_wcal,col_name,j,a);
00384        }
00385      }
00386      for (i=0; i< cfg->nrDispCoefficients; i++) {
00387        snprintf(col_name,MAX_NAME_SIZE-1,"%s%d","coeff",i);
00388        check_nomsg(coef_avg=cpl_table_get_column_mean(tbl_wcal,col_name));
00389        check_nomsg(coef_med=cpl_table_get_column_median(tbl_wcal,col_name));
00390 
00391        snprintf(key_name,MAX_NAME_SIZE-1,"%s%d%s","QC COEF",i," AVG");
00392        ck0_nomsg(sinfo_qclog_add_double(qclog_tbl,key_name,coef_avg,
00393                                         "Average wavecal Coef","%g"));
00394 
00395        snprintf(key_name,MAX_NAME_SIZE-1,"%s%d%s","QC COEF",i," MED");
00396        ck0_nomsg(sinfo_qclog_add_double(qclog_tbl,key_name,coef_med,
00397                                         "Median wavecal Coef","%g"));
00398 
00399      }
00400      if(pdensity >1) {
00401      strcpy(tbl_name,cfg->coeffsName);
00402      ck0(sinfo_pro_save_tbl(tbl_wcal,ref_set,sof,tbl_name,
00403                 PRO_WAVE_COEF_SLIT,qclog_tbl,plugin_id,config),
00404      "cannot save tbl %s", tbl_name);
00405      sinfo_free_table(&tbl_wcal);
00406      sinfo_free_table(&qclog_tbl);
00407      }
00408 
00409    }
00410 
00411    /*
00412     #store the resulting Gaussian fit parameters in an ASCII file if wished
00413    */
00414    if (cfg->writeParInd == 1) {
00415       sinfo_new_dump_fit_params_to_ascii(par,WAVECAL_FIT_PARAMS_OUT_FILEASCII);
00416  
00417       cknull(par,"no fit parameters available!") ;
00418       cknull(cfg->paramsList,"no filename available!") ;
00419       check_nomsg(tbl_fp = cpl_table_new(par[0] -> n_params));
00420       check_nomsg(cpl_table_new_column(tbl_fp,"n_params", CPL_TYPE_INT));
00421       check_nomsg(cpl_table_new_column(tbl_fp,"column", CPL_TYPE_INT));
00422       check_nomsg(cpl_table_new_column(tbl_fp,"line", CPL_TYPE_INT));
00423  
00424       for(j=0;j<4;j++) {
00425          snprintf(col,MAX_NAME_SIZE-1,"%s%d","fpar",j);
00426          cpl_table_new_column(tbl_fp,col, CPL_TYPE_DOUBLE);
00427          snprintf(col,MAX_NAME_SIZE-1,"%s%d","dpar",j);
00428          cpl_table_new_column(tbl_fp,col, CPL_TYPE_DOUBLE);
00429       }
00430 
00431       cknull_nomsg(qclog_tbl = sinfo_qclog_init());
00432       ck0_nomsg(sinfo_qclog_add_int(qclog_tbl,"QC NLINES",n_lines,
00433                                     "Number of found lines","%d"));
00434 
00435      for ( i = 0 ; i < par[0] -> n_params ; i++ ) {
00436       
00437          check_nomsg(cpl_table_set_int(tbl_fp,"n_params",i,par[i]->n_params));
00438          check_nomsg(cpl_table_set_int(tbl_fp,"column",i,par[i]->column));
00439          check_nomsg(cpl_table_set_int(tbl_fp,"line",i,par[i]->line));
00440 
00441          for(j=0;j<4;j++) {
00442        snprintf(col,MAX_NAME_SIZE-1,"%s%d","fpar",j);
00443        if(isnan(par[i]->fit_par[j])) {
00444          cpl_table_set_invalid(tbl_fp,col,i);
00445        } else {
00446          cpl_table_set_double(tbl_fp,col,i,par[i]->fit_par[j]);
00447        }
00448        snprintf(col,MAX_NAME_SIZE-1,"%s%d","dpar",j);
00449        if(isnan(par[i]->derv_par[j])) {
00450          cpl_table_set_invalid(tbl_fp,col,i);
00451        } else {
00452          cpl_table_set_double(tbl_fp,col,i,par[i]->derv_par[j]);
00453        }
00454      }
00455       }
00456 
00457       check_nomsg(fwhm_avg = cpl_table_get_column_mean(tbl_fp,"fpar1"));
00458       check_nomsg(fwhm_med = cpl_table_get_column_median(tbl_fp,"fpar1"));
00459       ck0_nomsg(sinfo_qclog_add_double(qclog_tbl,"QC FWHM MED",fwhm_med,
00460                                        "Median FWHM of found lines","%f"));
00461       ck0_nomsg(sinfo_qclog_add_double(qclog_tbl,"QC FWHM AVG",fwhm_avg,
00462                                        "Average FWHM of found lines","%f"));
00463 
00464       if(pdensity > 1) {
00465       ck0(sinfo_pro_save_tbl(tbl_fp,ref_set,sof,cfg->paramsList,
00466                  PRO_WAVE_PAR_LIST,qclog_tbl,plugin_id,config),
00467       "cannot save tbl %s", cfg->paramsList);
00468       }
00469 
00470       sinfo_free_table(&qclog_tbl);
00471       sinfo_free_table(&tbl_fp) ;
00472 
00473    }
00474 
00475    /* free memory */
00476    sinfo_new_destroy_2Dfloatarray ( &wavelength_clean, lx );
00477    sinfo_new_destroy_2Dintarray (&row_clean, lx);
00478    sinfo_new_destroy_intarray(&n_found_lines );
00479    sinfo_new_destroy_intarray(&sum_pointer );
00480    sinfo_new_destroy_2Dfloatarray ( &acoefs, cfg->nrDispCoefficients );
00481    sinfo_free_table(&tbl_line_list);
00482 
00483    /*----------------------------------------------------------------------
00484     *-------------------WAVEMAP--------------------------------------------
00485     *----------------------------------------------------------------------
00486     */
00487 
00488    /*
00489     #---now do the cross sinfo_correlation and produce a wavelength map---
00490    */
00491 
00492   } else if (cfg->wavemapInd == 1 && cfg->calibIndicator == 0) { 
00493     sinfo_msg("Wavemap");
00494     acoefs = sinfo_new_2Dfloatarray ( cfg->nrDispCoefficients, lx);
00495     /* #read the parameterized dispersion relation */
00496 
00497     strcpy(tbl_name,cfg->coeffsName);
00498     check_nomsg(tbl_wcal = cpl_table_load(tbl_name,1,0));
00499 
00500     for (i =0; i < lx; i++) {
00501       for (j = 0; j< cfg->nrDispCoefficients; j++) {
00502     snprintf(col_name,MAX_NAME_SIZE-1,"%s%d","coeff",j);
00503     acoefs[j][i]=cpl_table_get_double(tbl_wcal,col_name,i,status);
00504       }
00505     }
00506     sinfo_free_table(&tbl_wcal);
00507 
00508     cknull(map_img = sinfo_new_create_shifted_slit_wavemap2(im,
00509                            acoefs,
00510                            cfg->nrDispCoefficients,
00511                            wave,
00512                            intens,
00513                            n_lines,
00514                            cfg->magFactor,
00515                            cfg->guessDispersion1,
00516                            cfg->pixeldist ),
00517            "sinfo_createShiftedSlitWavemap2 failed!");
00518 
00519     par = sinfo_new_fit_params(15*n_lines);
00520     sinfo_msg("Check shifts");
00521 
00522     shift = sinfo_new_check_correlated_line_positions ( im, acoefs, 
00523                          cfg->nrDispCoefficients, 
00524                          wave,
00525                          intens,
00526                          n_lines,
00527                          cfg->fwhm,
00528                          cfg->halfWidth,
00529                          cfg->minAmplitude,
00530                          cfg->guessDispersion1,
00531                          par );
00532 
00533     if (FLAG == shift){
00534       sinfo_msg_error("sinfo_checkCorrelatedLinePositions failed!\n");
00535     }
00536 
00537     sinfo_det_ncounts(raw, cfg->qc_thresh_max,qc);
00538     cknull_nomsg(qclog_tbl = sinfo_qclog_init());
00539 
00540     ck0_nomsg(sinfo_qclog_add_int(qclog_tbl,"QC NLINES",n_lines,
00541                                   "Number of found lines","%d"));
00542 
00543     check_nomsg(fwhm_avg = cpl_table_get_column_mean(tbl_fp,"fpar1"));
00544     check_nomsg(fwhm_med = cpl_table_get_column_median(tbl_fp,"fpar1"));
00545 
00546     ck0_nomsg(sinfo_qclog_add_double(qclog_tbl,"QC FWHM MED",fwhm_med,
00547                                   "Median FWHM of found lines","%f"));
00548     ck0_nomsg(sinfo_qclog_add_double(qclog_tbl,"QC FWHM AVG",fwhm_avg,
00549                                   "Average FWHM of found lines","%f"));
00550 
00551     ck0(sinfo_pro_save_ima(map_img,ref_set,sof,cfg->outName,
00552                PRO_WAVE_MAP,qclog_tbl,plugin_id,config),
00553     "cannot save ima %s", cfg->outName);
00554 
00555     sinfo_free_table(&qclog_tbl);
00556     sinfo_free_image(&map_img);
00557 
00558     /* # ---free memory--- */
00559     sinfo_new_destroy_2Dfloatarray ( &acoefs, cfg->nrDispCoefficients );
00560     /* To fix a memory bug we comment the following. But check! */
00561     /* cpl_free ( wave ); */
00562     /* cpl_free ( intens ); */
00563 
00564   } else if (cfg->wavemapInd == 1 && cfg->calibIndicator == 1) {
00565     sinfo_msg_error("give either wavemapIndicator = yes and calibIndicator");
00566     sinfo_msg_error("= no or wavemapIndicator = no and calibIndicator = yes") ;
00567     goto cleanup;
00568   }
00569 
00570   /*-------------------------------------------------------------------
00571    *-------------------SLITFITS----------------------------------------
00572    *-------------------------------------------------------------------
00573    #--fit the slitlet sinfo_edge positions if desired--
00574   */
00575   if (cfg->slitposIndicator == 1) {
00576   sinfo_msg("fit the slitlet sinfo_edge positions");
00577   if (cfg->calibIndicator != 1 && cfg->estimateIndicator != 1) {
00578 
00579     /*
00580       #read the first integer to determine the number of fit 
00581       parameters to allocate
00582     */
00583 
00584      if(sinfo_is_fits_file(cfg->paramsList) !=1 ) {
00585        sinfo_msg_error("cannot read FITS file %s ",cfg->paramsList);
00586        goto cleanup;
00587      }
00588      strcpy(tbl_fitpar_name,cfg->paramsList);
00589      check_nomsg(tbl_fitpar = cpl_table_load(tbl_fitpar_name,1,0));
00590 
00591      check_nomsg(readsum=cpl_table_get_int(tbl_fitpar,"n_params",1,status));
00592      sinfo_free_table(&tbl_fitpar);
00593  
00594      cknull(sinfo_new_fit_params( readsum ),"sinfo_new_fit_params failed!");
00595 
00596      ck0(sinfo_dumpTblToFitParams(par, cfg->paramsList),
00597      "reading tbl %s ", cfg->paramsList);
00598   }
00599 
00600   /* #allocate memory for the slitlet position array */
00601   sinfo_slit_pos = sinfo_new_2Dfloatarray(32,2);
00602   /* #now decide which fit model function you want to use by 
00603      reading a given character
00604   */
00605   /*
00606   sinfo_msg("cfg->fitBoltzIndicator=%d\n",cfg->fitBoltzIndicator);
00607   sinfo_msg("cfg->estimateIndicator=%d\n",cfg->estimateIndicator);
00608   sinfo_msg("cfg->fitEdgeIndicator=%d\n",cfg->fitEdgeIndicator);
00609   */
00610   if (cfg->fitBoltzIndicator == 1) {
00611     if (cfg->estimateIndicator == 1) {
00612       /* #open the ASCII list of the slitlet positions--- */
00613       /*READ TFITS TABLE*/
00614       strcpy(tbl_slitpos_guess_name,cfg->slitposGuessName);
00615       check_nomsg(tbl_slitpos_guess=cpl_table_load(tbl_slitpos_guess_name,1,0));
00616       check_nomsg(n = cpl_table_get_nrow(tbl_slitpos_guess));
00617 
00618       for (i =0 ; i< 32; i++){
00619         val_x=cpl_table_get_double(tbl_slitpos_guess,"pos1",i,status);
00620         val_y=cpl_table_get_double(tbl_slitpos_guess,"pos2",i,status);
00621     sinfo_new_array2D_set_value(sinfo_slit_pos,val_x,i,0);
00622     sinfo_new_array2D_set_value(sinfo_slit_pos,val_y,i,1);
00623       }
00624       sinfo_free_table(&tbl_slitpos_guess);
00625 
00626       sinfo_msg("sinfo_new_fit_slits_boltz_with_estimate");
00627       fit = sinfo_new_fit_slits_boltz_with_estimate(im, 
00628                                             sinfo_slit_pos, 
00629                                   cfg->boxLength, 
00630                                             cfg->yBox, 
00631                                             cfg->diffTol, 
00632                                             cfg->loPos, 
00633                                             cfg->hiPos );
00634       if (fit < 0){
00635     sinfo_msg_error("sinfo_new_fit_slits_boltz_with_estimate failed" );
00636     goto cleanup;
00637       }
00638     } else {
00639     sinfo_msg("sinfo_new_fit_slits_boltz");
00640     fit = sinfo_new_fit_slits_boltz(im, 
00641                                   par, 
00642                                   sinfo_slit_pos, 
00643                                   cfg->boxLength, 
00644                                   cfg->yBox, 
00645                                   cfg->diffTol );
00646 
00647     if (fit < 0) {
00648       sinfo_msg_error ( "sinfo_new_fit_slits_boltz failed" );
00649       goto cleanup;
00650     }
00651     }
00652   } else if (cfg->fitEdgeIndicator == 1) {
00653 
00654     if (cfg->estimateIndicator == 1){
00655       /*READ TFITS TABLE*/
00656       strcpy(tbl_slitpos_guess_name,cfg->slitposGuessName);
00657       check_nomsg(tbl_slitpos_guess=cpl_table_load(tbl_slitpos_guess_name,1,0));
00658       check_nomsg(n = cpl_table_get_nrow(tbl_slitpos_guess));
00659 
00660       for (i =0 ; i< 32; i++){
00661     val_x=cpl_table_get_double(tbl_slitpos_guess,"pos1",i,status);
00662     val_y=cpl_table_get_double(tbl_slitpos_guess,"pos2",i,status);
00663     sinfo_new_array2D_set_value(sinfo_slit_pos,val_x,i,0);
00664     sinfo_new_array2D_set_value(sinfo_slit_pos,val_y,i,1);
00665       }
00666       cpl_table_delete(tbl_slitpos_guess);
00667 
00668       sinfo_msg("sinfo_new_fit_slits_edge_with_estimate");
00669       fit = sinfo_new_fit_slits_edge_with_estimate(im, 
00670                                            sinfo_slit_pos,
00671                                            cfg->boxLength, 
00672                                            cfg->yBox, 
00673                                            cfg->diffTol, 
00674                                            cfg->loPos, 
00675                                  cfg->hiPos );
00676       if (fit < 0) {
00677      sinfo_msg_error( "sinfo_new_fit_slits_boltz failed" );
00678          goto cleanup;
00679       }
00680     } else {
00681       sinfo_msg("sinfo_new_fit_slits_edge");
00682       fit = sinfo_new_fit_slits_edge(im, 
00683                                par, 
00684                                sinfo_slit_pos, 
00685                                cfg->boxLength, 
00686                                cfg->yBox, 
00687                                cfg->diffTol );
00688       if (fit < 0) {
00689     sinfo_msg_error("sinfo_new_fit_slits_edge failed" );
00690         goto cleanup;
00691       }
00692     }
00693   } else {
00694     sinfo_msg_error("no indication of desired fit function given" );
00695     goto cleanup;
00696   }
00697   sinfo_free_image(&im);
00698 
00699   /* #store the resulting sitlet positions in an TFITS table */
00700   check_nomsg(tbl_spos = cpl_table_new(32));
00701   check_nomsg(cpl_table_new_column(tbl_spos,"pos1", CPL_TYPE_DOUBLE));
00702   check_nomsg(cpl_table_new_column(tbl_spos,"pos2", CPL_TYPE_DOUBLE));
00703   check_nomsg(cpl_table_set_column_format(tbl_spos,"pos1", "15.9f"));
00704   check_nomsg(cpl_table_set_column_format(tbl_spos,"pos2", "15.9f"));
00705 
00706   for (i =0; i< 32; i++) {
00707     cpl_table_set_double(tbl_spos,"pos1",i,
00708                          sinfo_new_array2D_get_value(sinfo_slit_pos,i,0));
00709     cpl_table_set_double(tbl_spos,"pos2",i,
00710                          sinfo_new_array2D_get_value(sinfo_slit_pos,i,1));
00711      
00712   }
00713   if(sw == 1) {
00714     strcpy(tbl_name,"out_slitpos_guess.fits");
00715     ck0(sinfo_pro_save_tbl(tbl_spos,ref_set,sof,tbl_name,
00716                PRO_SLIT_POS_GUESS,NULL,plugin_id,config),
00717     "cannot save tbl %s", tbl_name);
00718   } else {
00719     strcpy(tbl_name,cfg->slitposName);
00720     ck0(sinfo_pro_save_tbl(tbl_spos,ref_set,sof,tbl_name,
00721                PRO_SLIT_POS,NULL,plugin_id,config),
00722     "cannot save tbl %s", tbl_name);
00723   }
00724   sinfo_free_table(&tbl_spos);
00725   sinfo_new_destroy_2Dfloatarray ( &sinfo_slit_pos, 32 );
00726   }
00727 
00728   if ( (cfg->slitposIndicator == 1 && cfg->estimateIndicator != 1) || 
00729        (cfg->calibIndicator == 1)  || (cfg->wavemapInd == 1) ){
00730     sinfo_new_destroy_fit_params(&par);
00731   }
00732 
00733 
00734 
00735 
00736 
00737   if(pdensity > 1) {
00738      check_nomsg(sinfo_image_resample(plugin_id,config,sof,ref_set));
00739 
00740 
00796   }
00797   sinfo_free_frameset(&raw);
00798   sinfo_qc_wcal_delete(&qc);
00799   sinfo_wavecal_free(&cfg);
00800 
00801   return 0;
00802 
00803   cleanup:
00804   sinfo_free_image(&map_img);
00805   //sinfo_free_image(&wstk_img);
00806   sinfo_free_table(&tbl_spos);
00807   sinfo_free_table(&tbl_fitpar);
00808   sinfo_free_image(&map_img);
00809   sinfo_free_table(&tbl_wcal);
00810   sinfo_free_table(&tbl_fp) ;
00811   sinfo_free_table(&tbl_line_list);
00812   sinfo_free_table(&tbl_wcal);
00813   sinfo_free_table(&qclog_tbl);
00814   sinfo_free_image(&map_img);
00815   sinfo_new_destroy_fit_params(&par);
00816   sinfo_new_destroy_2Dfloatarray ( &sinfo_slit_pos, 32 );
00817   sinfo_new_destroy_2Dfloatarray ( &wavelength_clean, lx );
00818   sinfo_new_destroy_2Dintarray (&row_clean, lx);
00819   sinfo_new_destroy_intarray(&n_found_lines );
00820   sinfo_new_destroy_intarray(&sum_pointer );
00821   if(acoefs!=NULL) {
00822      sinfo_new_destroy_2Dfloatarray(&acoefs,cfg->nrDispCoefficients);
00823   }
00824   sinfo_free_table(&tbl_line_list);
00825   sinfo_free_image(&im);
00826   sinfo_wavecal_free(&cfg);
00827   sinfo_free_frameset(&raw);
00828   sinfo_qc_wcal_delete(&qc);
00829   return -1 ;
00830 
00831 }
00832 
00833 
00834 
00835 
00836 
00837 
00838 
00839 
00840 
00841 
00842 
00843 
00844 static cpl_error_code
00845 sinfo_image_resample(const char* plugin_id,
00846                      cpl_parameterlist* config,
00847                      cpl_frameset* sof,
00848                      cpl_frameset* ref_set)
00849 {
00850   //RESAMPLE ThAr frame for QC
00851   double dis=0;
00852   float mi=0;
00853   float ma=0;
00854   double cwav=0;
00855   int cpix=0;
00856   const cpl_frame* frm=NULL;
00857   char wstk_name[80];
00858   char map_name[80];
00859   cpl_image* res_ima=NULL;
00860   int ncoeffs=3;
00861   int nrows=SINFO_RESAMP_NROWS;
00862   cpl_parameter* p=NULL;
00863   cpl_image* wstk_img=NULL;
00864   cpl_image* map_img=NULL;
00865 
00866   check_nomsg(p = cpl_parameterlist_find(config,
00867                                          "sinfoni.wavecal.n_coeffs"));
00868 
00869   check_nomsg(ncoeffs=cpl_parameter_get_int(p));
00870 
00871   check_nomsg(frm=cpl_frameset_find_const(sof,PRO_WAVE_LAMP_STACKED));
00872   check_nomsg(strcpy(wstk_name,cpl_frame_get_filename(frm)));
00873   check_nomsg(wstk_img=cpl_image_load(wstk_name,CPL_TYPE_FLOAT,0,0));
00874 
00875 
00876 
00877   check_nomsg(frm=cpl_frameset_find_const(sof,PRO_WAVE_MAP));
00878   check_nomsg(strcpy(map_name,cpl_frame_get_filename(frm)));
00879   check_nomsg(map_img=cpl_image_load(map_name,CPL_TYPE_FLOAT,0,0));
00880 
00881 
00882 
00883   cknull(res_ima = sinfo_new_defined_resampling(wstk_img, 
00884                         map_img, 
00885                         ncoeffs,
00886                         &nrows,
00887                         &dis,
00888                         &mi,
00889                         &ma,
00890                         &cwav,
00891                         &cpix),
00892      " sinfo_definedResampling() failed" ) ;
00893 
00894 
00895 
00896   ck0(sinfo_pro_save_ima(res_ima,ref_set,sof,WAVECAL_RESAMPLED_OUT_FILENAME,
00897              PRO_RESAMPLED_WAVE,NULL,plugin_id,config),
00898       "cannot save ima %s",WAVECAL_RESAMPLED_OUT_FILENAME);
00899 
00900 
00901   cleanup:
00902 
00903   sinfo_free_image(&map_img);
00904   sinfo_free_image(&res_ima);
00905   sinfo_free_image(&wstk_img);
00906   return cpl_error_get_code();
00907 
00908 
00909 }

Generated on 3 Mar 2013 for SINFONI Pipeline Reference Manual by  doxygen 1.6.1