sinfo_new_resampling.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    
00022    File name     :    sinfo_new_resampling.c
00023    Author         :    Nicolas Devillard
00024    Created on    :    Jan 04, 1996
00025    Description    :    resampling routines
00026 
00027  ---------------------------------------------------------------------------*/
00028 
00029 /*
00030     $Id: sinfo_new_resampling.c,v 1.9 2007/06/06 07:10:45 amodigli Exp $
00031     $Author: amodigli $
00032     $Date: 2007/06/06 07:10:45 $
00033     $Revision: 1.9 $
00034  */
00035 
00036 #ifdef HAVE_CONFIG_H
00037 #  include <config.h>
00038 #endif
00039 /*---------------------------------------------------------------------------
00040                                   Includes
00041  ---------------------------------------------------------------------------*/
00042 #include <math.h>
00043 #include "sinfo_new_resampling.h"
00044 #include "sinfo_pixel_handling.h"
00045 #include "sinfo_globals.h"
00046 /* #include "my_pi.h" */
00047 #include "sinfo_resampling.h"
00048 /*---------------------------------------------------------------------------
00049                               Private functions
00050  ---------------------------------------------------------------------------*/
00051 
00052 static void new_reverse_tanh_kernel(double * data, int nn) ;
00053 static double sinfo_new_sinc(double x);
00054 
00055 /*---------------------------------------------------------------------------
00056                               Function codes
00057  ---------------------------------------------------------------------------*/
00058 
00066 /*-------------------------------------------------------------------------*/
00090 /*--------------------------------------------------------------------------*/
00091 
00092 double   *
00093 sinfo_new_generate_interpolation_kernel(const char * kernel_type)
00094 {
00095     double  *    tab ;
00096     int         i ;
00097     double      x ;
00098     double        alpha ;
00099     double        inv_norm ;
00100     int         samples = KERNEL_SAMPLES ;
00101 
00102     if (kernel_type==NULL) {
00103         tab = sinfo_new_generate_interpolation_kernel("tanh") ;
00104     } else if (!strcmp(kernel_type, "default")) {
00105         tab = sinfo_new_generate_interpolation_kernel("tanh") ;
00106     } else if (!strcmp(kernel_type, "sinfo_new_sinc")) {
00107         tab = cpl_malloc(samples * sizeof(double)) ;
00108         tab[0] = 1.0 ;
00109         tab[samples-1] = 0.0 ;
00110         for (i=1 ; i<samples ; i++) {
00111             x = (double)KERNEL_WIDTH * (double)i/(double)(samples-1) ;
00112             tab[i] = sinfo_new_sinc(x) ;
00113         }
00114     } else if (!strcmp(kernel_type, "sinc2")) {
00115         tab = cpl_malloc(samples * sizeof(double)) ;
00116         tab[0] = 1.0 ;
00117         tab[samples-1] = 0.0 ;
00118         for (i=1 ; i<samples ; i++) {
00119             x = 2.0 * (double)i/(double)(samples-1) ;
00120             tab[i] = sinfo_new_sinc(x) ;
00121             tab[i] *= tab[i] ;
00122         }
00123     } else if (!strcmp(kernel_type, "lanczos")) {
00124         tab = cpl_malloc(samples * sizeof(double)) ;
00125         for (i=0 ; i<samples ; i++) {
00126             x = (double)KERNEL_WIDTH * (double)i/(double)(samples-1) ;
00127             if (fabs(x)<2) {
00128                 tab[i] = sinfo_new_sinc(x) * sinfo_new_sinc(x/2) ;
00129             } else {
00130                 tab[i] = 0.00 ;
00131             }
00132         }
00133     } else if (!strcmp(kernel_type, "hamming")) {
00134         tab = cpl_malloc(samples * sizeof(double)) ;
00135         alpha = 0.54 ;
00136         inv_norm  = 1.00 / (double)(samples - 1) ;
00137         for (i=0 ; i<samples ; i++) {
00138             x = (double)i ;
00139             if (i<(samples-1)/2) {
00140                 tab[i] = alpha + (1-alpha) * cos(2.0*PI_NUMB*x*inv_norm) ;
00141             } else {
00142                 tab[i] = 0.0 ;
00143             }
00144         }
00145     } else if (!strcmp(kernel_type, "hann")) {
00146         tab = cpl_malloc(samples * sizeof(double)) ;
00147         alpha = 0.50 ;
00148         inv_norm  = 1.00 / (double)(samples - 1) ;
00149         for (i=0 ; i<samples ; i++) {
00150             x = (double)i ;
00151             if (i<(samples-1)/2) {
00152                 tab[i] = alpha + (1-alpha) * cos(2.0*PI_NUMB*x*inv_norm) ;
00153             } else {
00154                 tab[i] = 0.0 ;
00155             }
00156         }
00157     } else if (!strcmp(kernel_type, "tanh")) {
00158         tab = sinfo_new_generate_tanh_kernel(TANH_STEEPNESS) ;
00159     } else {
00160         sinfo_msg_error("unrecognized kernel type [%s]: aborting generation",
00161                 kernel_type) ;
00162         return NULL ;
00163     }
00164 
00165     return tab ;
00166 }
00167 
00168 /*-------------------------------------------------------------------------*/
00179 /*--------------------------------------------------------------------------*/
00180 
00181 double
00182 sinfo_new_sinc(double x)
00183 {
00184     if (fabs(x)<1e-4)
00185         return (double)1.00 ;
00186     else
00187         return ((sin(x * (double)PI_NUMB)) / (x * (double)PI_NUMB)) ;
00188 }
00189 
00190 
00191 
00192 /*-------------------------------------------------------------------------*/
00228 /*--------------------------------------------------------------------------*/
00229 
00230 cpl_image *
00231 sinfo_new_warp_image_generic(
00232     cpl_image         *    image_in,
00233     char          *    kernel_type,
00234     cpl_polynomial    *    poly_u,
00235     cpl_polynomial    *    poly_v
00236 )
00237 {
00238     cpl_image    *    image_out ;
00239     int             i, j, k ;
00240     int             lx_out, ly_out ;
00241     double           cur ;
00242     double           neighbors[16] ;
00243     double           rsc[8],
00244                     sumrs ;
00245     double           x, y ;
00246     int             px, py ;
00247     int             pos ;
00248     int             tabx, taby ;
00249     double      *    kernel ;
00250     int                  leaps[16] ;
00251     int ilx=0;
00252     int ily=0;
00253     float* pidata=NULL;
00254     float* podata=NULL;
00255     cpl_vector* vx=NULL;
00256     if (image_in == NULL) return NULL ;
00257 
00258 
00259     /* Generate default interpolation kernel */
00260     kernel = sinfo_new_generate_interpolation_kernel(kernel_type) ;
00261     if (kernel == NULL) {
00262         sinfo_msg_error("cannot generate kernel: aborting resampling") ;
00263         return NULL ;
00264     }
00265     ilx=cpl_image_get_size_x(image_in);
00266     ily=cpl_image_get_size_y(image_in);
00267     pidata=cpl_image_get_data_float(image_in);
00268 
00269     /* Compute new image size   */
00270     lx_out = (int)ilx ;
00271     ly_out = (int)ily ;
00272 
00273     image_out = cpl_image_new(lx_out, ly_out,CPL_TYPE_FLOAT) ;
00274     podata=cpl_image_get_data_float(image_out);
00275 
00276     /* Pre compute leaps for 16 closest neighbors positions */
00277 
00278     leaps[0] = -1 - ilx ;
00279     leaps[1] =    - ilx ;
00280     leaps[2] =  1 - ilx ;
00281     leaps[3] =  2 - ilx ;
00282 
00283     leaps[4] = -1 ;
00284     leaps[5] =  0 ;
00285     leaps[6] =  1 ;
00286     leaps[7] =  2 ;
00287 
00288     leaps[8] = -1 + ilx ;
00289     leaps[9] =      ilx ;
00290     leaps[10]=  1 + ilx ;
00291     leaps[11]=  2 + ilx ;
00292 
00293     leaps[12]= -1 + 2*ilx ;
00294     leaps[13]=      2*ilx ;
00295     leaps[14]=  1 + 2*ilx ;
00296     leaps[15]=  2 + 2*ilx ;
00297 
00298     vx=cpl_vector_new(2); /* vector of size 2 = polynomial dimension */
00299     /* Double loop on the output image  */
00300     for (j=0 ; j < ly_out ; j++) {
00301         for (i=0 ; i< lx_out ; i++) {
00302             /* Compute the original source for this pixel   */
00303       cpl_vector_set(vx,0,(double)i);
00304       cpl_vector_set(vx,1,(double)j);
00305             x = cpl_polynomial_eval(poly_u, vx);
00306             y = cpl_polynomial_eval(poly_v, vx);
00307 
00308             /* Which is the closest integer positioned neighbor?    */
00309             px = (int)x ;
00310             py = (int)y ;
00311 
00312             if ((px < 1) ||
00313                 (px > (ilx-3)) ||
00314                 (py < 1) ||
00315                 (py > (ily-3)))
00316                 podata[i+j*lx_out] = (pixelvalue)0.0/0.0 ;
00317             else {
00318                 /* Now feed the positions for the closest 16 neighbors  */
00319                 pos = px + py * ilx ;
00320                 for (k=0 ; k<16 ; k++)
00321                     neighbors[k] = (double)(pidata[(int)(pos+leaps[k])]) ;
00322 
00323                 /* Which tabulated value index shall we use?    */
00324                 tabx = (x - (double)px) * (double)(TABSPERPIX) ; 
00325                 taby = (y - (double)py) * (double)(TABSPERPIX) ; 
00326 
00327                 /* Compute resampling coefficients  */
00328                 /* rsc[0..3] in x, rsc[4..7] in y   */
00329 
00330                 rsc[0] = kernel[TABSPERPIX + tabx] ;
00331                 rsc[1] = kernel[tabx] ;
00332                 rsc[2] = kernel[TABSPERPIX - tabx] ;
00333                 rsc[3] = kernel[2 * TABSPERPIX - tabx] ;
00334                 rsc[4] = kernel[TABSPERPIX + taby] ;
00335                 rsc[5] = kernel[taby] ;
00336                 rsc[6] = kernel[TABSPERPIX - taby] ;
00337                 rsc[7] = kernel[2 * TABSPERPIX - taby] ;
00338 
00339                 sumrs = (rsc[0]+rsc[1]+rsc[2]+rsc[3]) *
00340                         (rsc[4]+rsc[5]+rsc[6]+rsc[7]) ;
00341 
00342                 /* Compute interpolated pixel now   */
00343                 cur =   rsc[4] * (  rsc[0]*neighbors[0] +
00344                                     rsc[1]*neighbors[1] +
00345                                     rsc[2]*neighbors[2] +
00346                                     rsc[3]*neighbors[3] ) +
00347                         rsc[5] * (  rsc[0]*neighbors[4] +
00348                                     rsc[1]*neighbors[5] +
00349                                     rsc[2]*neighbors[6] +
00350                                     rsc[3]*neighbors[7] ) +
00351                         rsc[6] * (  rsc[0]*neighbors[8] +
00352                                     rsc[1]*neighbors[9] +
00353                                     rsc[2]*neighbors[10] +
00354                                     rsc[3]*neighbors[11] ) +
00355                         rsc[7] * (  rsc[0]*neighbors[12] +
00356                                     rsc[1]*neighbors[13] +
00357                                     rsc[2]*neighbors[14] +
00358                                     rsc[3]*neighbors[15] ) ; 
00359 
00360                 /* Affect the value to the output image */
00361                 podata[i+j*lx_out] = (pixelvalue)(cur/sumrs) ;
00362                 /* done ! */
00363             }       
00364         }
00365     }
00366     cpl_vector_delete(vx);
00367     cpl_free(kernel) ;
00368     return image_out ;
00369 }
00370 
00371 
00372 /*-------------------------------------------------------------------------*/
00393 /*--------------------------------------------------------------------------*/
00394 
00395 #define hk_gen(x,s) (((tanh(s*(x+0.5))+1)/2)*((tanh(s*(-x+0.5))+1)/2))
00396 
00397 double * sinfo_new_generate_tanh_kernel(double steep)
00398 {
00399     double  *   kernel ;
00400     double  *   x ;
00401     double      width ;
00402     double      inv_np ;
00403     double      ind ;
00404     int         i ;
00405     int         np ;
00406     int         samples ;
00407 
00408     width   = (double)TABSPERPIX / 2.0 ; 
00409     samples = KERNEL_SAMPLES ;
00410     np      = 32768 ; /* Hardcoded: should never be changed */
00411     inv_np  = 1.00 / (double)np ;
00412 
00413     /*
00414      * Generate the kernel expression in Fourier space
00415      * with a correct frequency ordering to allow standard FT
00416      */
00417     x = cpl_malloc((2*np+1)*sizeof(double)) ;
00418     for (i=0 ; i<np/2 ; i++) {
00419         ind      = (double)i * 2.0 * width * inv_np ;
00420         x[2*i]   = hk_gen(ind, steep) ;
00421         x[2*i+1] = 0.00 ;
00422     }
00423     for (i=np/2 ; i<np ; i++) {
00424         ind      = (double)(i-np) * 2.0 * width * inv_np ;
00425         x[2*i]   = hk_gen(ind, steep) ;
00426         x[2*i+1] = 0.00 ;
00427     }
00428 
00429     /* 
00430      * Reverse Fourier to come back to image space
00431      */
00432     new_reverse_tanh_kernel(x, np) ;
00433 
00434     /*
00435      * Allocate and fill in returned array
00436      */
00437     kernel = cpl_malloc(samples * sizeof(double)) ;
00438     for (i=0 ; i<samples ; i++) {
00439         kernel[i] = 2.0 * width * x[2*i] * inv_np ;
00440     }
00441     cpl_free(x) ;
00442     return kernel ;
00443 }
00444 
00445 
00446 /*-------------------------------------------------------------------------*/
00458 /*--------------------------------------------------------------------------*/
00459 
00460 #define KERNEL_SW(a,b) tempr=(a);(a)=(b);(b)=tempr
00461 static void new_reverse_tanh_kernel(double * data, int nn)
00462 {
00463     unsigned long   n,
00464                     mmax,
00465                     m,
00466                     i, j,
00467                     istep ;
00468     double  wtemp,
00469             wr,
00470             wpr,
00471             wpi,
00472             wi,
00473             theta;
00474     double  tempr,
00475             tempi;
00476 
00477     n = (unsigned long)nn << 1;
00478     j = 1;
00479     for (i=1 ; i<n ; i+=2) {
00480         if (j > i) {
00481             KERNEL_SW(data[j-1],data[i-1]);
00482             KERNEL_SW(data[j],data[i]);
00483         }
00484         m = n >> 1;
00485         while (m>=2 && j>m) {
00486             j -= m;
00487             m >>= 1;
00488         }
00489         j += m;
00490     }
00491     mmax = 2;
00492     while (n > mmax) {
00493         istep = mmax << 1;
00494         theta = 2 * PI_NUMB / mmax;
00495         wtemp = sin(0.5 * theta);
00496         wpr = -2.0 * wtemp * wtemp;
00497         wpi = sin(theta);
00498         wr  = 1.0;
00499         wi  = 0.0;
00500         for (m=1 ; m<mmax ; m+=2) {
00501             for (i=m ; i<=n ; i+=istep) {
00502                 j = i + mmax;
00503                 tempr = wr * data[j-1] - wi * data[j];
00504                 tempi = wr * data[j]   + wi * data[j-1];
00505                 data[j-1] = data[i-1] - tempr;
00506                 data[j]   = data[i]   - tempi;
00507                 data[i-1] += tempr;
00508                 data[i]   += tempi;
00509             }
00510             wr = (wtemp = wr) * wpr - wi * wpi + wr;
00511             wi = wi * wpr + wtemp * wpi + wi;
00512         }
00513         mmax = istep;
00514     }
00515 }
00516 #undef KERNEL_SW
00517 
00518 
00519 
00520 /*-------------------------------------------------------------------------*/
00533 /*--------------------------------------------------------------------------*/
00534 
00535 void sinfo_new_show_interpolation_kernel(char * kernel_name)
00536 {
00537     double    *    ker ;
00538     int            i ;
00539     double        x ;
00540 
00541 
00542     ker = sinfo_new_generate_interpolation_kernel(kernel_name) ;
00543     if (ker == NULL)
00544         return ;
00545 
00546     (void)fprintf(stdout, "# Kernel is %s\n", kernel_name) ;
00547     x = 0.00 ;
00548     for (i=0 ; i<KERNEL_SAMPLES ; i++) {
00549         (void)fprintf(stdout, "%g %g\n", x, ker[i]) ;
00550         x += 1.00 / (double)TABSPERPIX ;
00551     }
00552     cpl_free(ker) ;
00553     return ;
00554 }

Generated on 8 Mar 2011 for SINFONI Pipeline Reference Manual by  doxygen 1.6.1