irplib_fft_body.h

00001 /* $Id: irplib_fft_body.h,v 1.1 2008/06/03 15:56:32 lbilbao Exp $
00002  *
00003  * This file is part of the irplib package 
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., 51 Franklin St, Fifth Floor, Boston, MA  02111-1307  USA
00019  */
00020 
00021 /*
00022  * $Author: lbilbao $
00023  * $Date: 2008/06/03 15:56:32 $
00024  * $Revision: 1.1 $
00025  * $Name: naco-4_1_2 $
00026  */
00027 
00028 #define CONCAT(a,b) a ## _ ## b
00029 #define CONCAT2X(a,b) CONCAT(a,b)
00030 
00031 #define IN_TYPE_ADD(a) CONCAT2X(a, IN_TYPE)
00032 #define IN_TYPE_CONST_ADD(a) CONCAT2X(IN_TYPE_ADD(a), const)
00033 
00034 #define OUT_TYPE_ADD(a) CONCAT2X(a, OUT_TYPE)
00035 
00036 #define TYPES_ADD(a) CONCAT2X(a, CONCAT2X(IN_TYPE, OUT_TYPE))
00037 
00038 #define CPL_TYPE_ADDED CONCAT2X(CPL_TYPE, FFTW_PREC_STRING)
00039 
00040 cpl_error_code TYPES_ADD(irplib_imagelist_fft)(cpl_imagelist       * real_out,
00041                            cpl_imagelist       * imag_out,
00042                            const cpl_imagelist * real_in,
00043                            const cpl_imagelist * imag_in,
00044                            unsigned              mode)
00045 {
00046     const cpl_image   * ref;
00047 
00048     fftw_complex      * in;
00049     fftw_complex      * out;
00050 
00051     unsigned            fftwmode;
00052 
00053     int                 nx, ny, nimg;
00054     int                 i, j;
00055 
00056     if (mode == IRPLIB_FFT_DEFAULT) {
00057     fftwmode = FFTW_FORWARD;
00058     } else if (mode == IRPLIB_FFT_INVERSE) {
00059     fftwmode = FFTW_BACKWARD;
00060     }
00061 
00062     ref  = cpl_imagelist_get_const(real_in, 0);
00063 
00064     nx   = cpl_image_get_size_x(ref);
00065     ny   = cpl_image_get_size_y(ref);
00066 
00067     nimg = cpl_imagelist_get_size(real_in);
00068 
00069     in = (fftw_complex *) fftw_malloc(sizeof(fftw_complex) *
00070                       nx * ny * nimg);
00071     out = in;
00072 
00073     /*
00074      * We can find three sub-cases:
00075      *     1.- Real input  (imag_in  == NULL)
00076      *     2.- Real output (imag_out == NULL)
00077      *     3.- None of the above (both input and output complex). Includes
00078      *         in-place transform (real_in == real_out && imag_in == imag_out)
00079      */
00080 
00081     /* Data preparation for sub-cases 2 and 3 */
00082     if (imag_in != NULL && imag_out != NULL) {
00083 
00084     for(i = 0; i < nimg; i++) {
00085         const IN_TYPE * real_data = 
00086         IN_TYPE_CONST_ADD(cpl_image_get_data)
00087         (cpl_imagelist_get_const(real_in, i));
00088         
00089         const IN_TYPE * imag_data =
00090         IN_TYPE_CONST_ADD(cpl_image_get_data)
00091         (cpl_imagelist_get_const(imag_in, i));
00092         
00093         for(j = 0; j < nx * ny; j++) {
00094         *((FFTW_PREC complex *)in + i * nx * ny + j) = 
00095             ((FFTW_PREC)(*(real_data + j)) + 
00096              (FFTW_PREC)(*(imag_data + j)) * I);
00097         }
00098     }
00099     } else if (imag_out == NULL) {
00100 
00101     for(i = 0; i < nimg; i++) {
00102         const IN_TYPE * real_data = 
00103         IN_TYPE_CONST_ADD(cpl_image_get_data)
00104         (cpl_imagelist_get_const(real_in, i));
00105         
00106         const IN_TYPE * imag_data =
00107         IN_TYPE_CONST_ADD(cpl_image_get_data)
00108         (cpl_imagelist_get_const(imag_in, i));
00109         
00110         for(j = 0; j < nx * (ny / 2 + 1); j++) {
00111         *((FFTW_PREC complex *)in + i * nx * ny + j) = 
00112             ((FFTW_PREC)(*(real_data + j)) + 
00113              (FFTW_PREC)(*(imag_data + j)) * I);
00114         }
00115     }
00116     }
00117 
00118     /* Now we proceed with the specific processing of each sub-case */
00119 
00120     if (imag_in == NULL) {
00121     /* Sub-case 1 */
00122 
00123     rfftwnd_plan rp = 
00124         rfftw2d_create_plan(nx, ny, FFTW_FORWARD, FFTW_ESTIMATE);
00125         
00126     for(i = 0; i < nimg; i++) {
00127 
00128         /*
00129          * If the input pixel type is coherent with the precision
00130          * of the FFTW library available, the pixel buffer can be 
00131          * directly handled by FFTW. Otherwise, it must be casted.
00132          */
00133 
00134 #if FFTW_PREC_CODE == IN_TYPE_CODE
00135         const IN_TYPE * real_data = 
00136         IN_TYPE_CONST_ADD(cpl_image_get_data)
00137         (cpl_imagelist_get_const(real_in, i));
00138         
00139         rfftwnd_real_to_complex(rp, 1, 
00140                     (fftw_real *)real_data, 1, nx * ny,
00141                     out + i * nx * ny,      1, nx * ny);
00142 #else
00143         const cpl_image * fprec_img = 
00144         cpl_image_cast(cpl_imagelist_get_const(real_in, i),
00145                    CPL_TYPE_ADDED);
00146 
00147         const FFTW_PREC * real_data = 
00148         (FFTW_PREC *)cpl_image_get_data_const(fprec_img);
00149         
00150         rfftwnd_real_to_complex(rp, 1, 
00151                     (fftw_real *)real_data, 1, nx * ny,
00152                     out + i * nx * ny,      1, nx * ny);
00153         cpl_image_delete((cpl_image *)fprec_img);
00154 #endif
00155     }
00156     
00157     rfftwnd_destroy_plan(rp);    
00158     } else if (imag_out == NULL) {
00159     /* Sub-case 2 */
00160         
00161     rfftwnd_plan rp = 
00162         rfftw2d_create_plan(nx, ny, FFTW_BACKWARD, FFTW_ESTIMATE);
00163         
00164 
00165     for(i = 0; i < nimg; i++) {
00166 
00167         OUT_TYPE * real_outdata = 
00168         OUT_TYPE_ADD(cpl_image_get_data)
00169         (cpl_imagelist_get(real_out, i));
00170 
00171         /*
00172          * Similarly to sub-case 2, the output pixel type can only
00173          * by directly written by FFTW if it is coherent with the
00174          * precision of the FFTW library available. Otherwise a new
00175          * pixel buffer of the appropriate precision is created and
00176          * copied to the real output buffer after computation.
00177          */
00178 
00179 #if FFTW_PREC_CODE == OUT_TYPE_CODE
00180         rfftwnd_complex_to_real(rp, 1,
00181                     in + i * nx * ny,          1, nx * ny,
00182                     (fftw_real *)real_outdata, 1, nx * ny);
00183 #else
00184         cpl_image * fprec_img = 
00185         cpl_image_new(nx, ny, CPL_TYPE_ADDED);
00186         FFTW_PREC * fp_outdata = 
00187         (FFTW_PREC *)cpl_image_get_data(fprec_img);
00188                   
00189         rfftwnd_complex_to_real(rp, 1, 
00190                     in + i * nx * ny,        1, nx * ny,
00191                     (fftw_real *)fp_outdata, 1, nx * ny);
00192 
00193         for(j = 0; j < nx * ny; j++) {
00194         *(real_outdata + j) = (OUT_TYPE)(*(fp_outdata + j));
00195         }
00196 
00197         cpl_image_delete(fprec_img);
00198 #endif
00199     }
00200 
00201     rfftwnd_destroy_plan(rp);
00202     } else {
00203     /* Sub-case 3 */
00204 
00205     fftwnd_plan p =
00206         fftw2d_create_plan_specific(nx, ny, fftwmode, FFTW_ESTIMATE |
00207                     FFTW_IN_PLACE,
00208                     in, 1, out, 1);
00209     fftwnd(p, nimg, in, 1, nx * ny, out, 1, nx * ny);
00210     fftwnd_destroy_plan(p);
00211     }
00212 
00213     /* Output preparation for sub-cases 1 and 3 */
00214     if (imag_out != NULL && imag_in != NULL) {
00215 
00216     for(i = 0; i < nimg; i++) {
00217         
00218         OUT_TYPE * real_outdata = 
00219         OUT_TYPE_ADD(cpl_image_get_data)
00220         (cpl_imagelist_get(real_out, i));
00221 
00222         OUT_TYPE * imag_outdata =
00223         OUT_TYPE_ADD(cpl_image_get_data)
00224         (cpl_imagelist_get(imag_out, i));
00225 
00226         for(j = 0; j < nx * ny; j++) {
00227         *(real_outdata + j) = 
00228             (OUT_TYPE)((out + i * nx * ny + j) -> re);
00229 
00230         *(imag_outdata + j) = 
00231             (OUT_TYPE)((out + i * nx * ny + j) -> im);
00232         }
00233     }
00234     } else if (imag_in == NULL) {
00235 
00236     for(i = 0; i < nimg; i++) {
00237         
00238         OUT_TYPE * real_outdata = 
00239         OUT_TYPE_ADD(cpl_image_get_data)
00240         (cpl_imagelist_get(real_out, i));
00241 
00242         OUT_TYPE * imag_outdata =
00243         OUT_TYPE_ADD(cpl_image_get_data)
00244         (cpl_imagelist_get(imag_out, i));
00245 
00246         int k;
00247 
00248         for(j = 0; j < nx; j++) {
00249         *(real_outdata + j) = 
00250             (OUT_TYPE)((out + i * nx * ny + j) -> re);
00251 
00252         *(imag_outdata + j) = 
00253             (OUT_TYPE)((out + i * nx * ny + j) -> im);
00254         }
00255 
00256         for(j = 1; j < (ny / 2 + 1); j++) {
00257         for(k = 0; k < nx; k++) {
00258             *(real_outdata + j * nx + k) = 
00259             (OUT_TYPE)((out + i * nx * ny + j * nx + k) -> re);
00260 
00261             *(real_outdata + (ny - j) * nx + k) = 
00262             (OUT_TYPE)((out + i * nx * ny + j * nx + k) -> re);
00263 
00264             *(imag_outdata + j * nx + k) = 
00265             (OUT_TYPE)((out + i * nx * ny + j * nx + k) -> im);
00266 
00267             *(imag_outdata + (ny - j) * nx + k) = 
00268             (OUT_TYPE)((out + i * nx * ny + j * nx + k) -> im);
00269         }
00270         }
00271     }
00272     }
00273 
00274     fftw_free(in);
00275     
00276     return cpl_error_get_code();
00277 }
00278 
00279 cpl_error_code TYPES_ADD(irplib_image_fft)(cpl_image       * real_out,
00280                        cpl_image       * imag_out,
00281                        const cpl_image * real_in,
00282                        const cpl_image * imag_in,
00283                        unsigned          mode)
00284 {
00285     fftw_complex      * in;
00286     fftw_complex      * out;
00287 
00288     unsigned            fftwmode;
00289 
00290     int                 nx, ny;
00291     int                 j;
00292 
00293     if (mode == IRPLIB_FFT_DEFAULT) {
00294     fftwmode = FFTW_FORWARD;
00295     } else if (mode == IRPLIB_FFT_INVERSE) {
00296     fftwmode = FFTW_BACKWARD;
00297     }
00298 
00299     nx   = cpl_image_get_size_x(real_in);
00300     ny   = cpl_image_get_size_y(real_in);
00301 
00302     in = (fftw_complex *) fftw_malloc(sizeof(fftw_complex) * nx * ny);
00303     out = in;
00304 
00305     /*
00306      * We can find three sub-cases:
00307      *     1.- Real input  (imag_in  == NULL)
00308      *     2.- Real output (imag_out == NULL)
00309      *     3.- None of the above (both input and output complex). Includes
00310      *         in-place transform (real_in == real_out && imag_in == imag_out)
00311      */
00312 
00313     /* Data preparation for sub-cases 2 and 3 */
00314     if (imag_in != NULL && imag_out != NULL) {
00315 
00316     const IN_TYPE * real_data = 
00317         IN_TYPE_CONST_ADD(cpl_image_get_data)(real_in);
00318         
00319     const IN_TYPE * imag_data =
00320         IN_TYPE_CONST_ADD(cpl_image_get_data)(imag_in);
00321         
00322     for(j = 0; j < nx * ny; j++) {
00323         *((FFTW_PREC complex *)in + j) = 
00324         ((FFTW_PREC)(*(real_data + j)) + 
00325          (FFTW_PREC)(*(imag_data + j)) * I);
00326     }
00327     } else if (imag_out == NULL) {
00328 
00329     const IN_TYPE * real_data = 
00330         IN_TYPE_CONST_ADD(cpl_image_get_data)(real_in);
00331         
00332     const IN_TYPE * imag_data =
00333         IN_TYPE_CONST_ADD(cpl_image_get_data)(imag_in);
00334         
00335     for(j = 0; j < nx * (ny / 2 + 1); j++) {
00336         *((FFTW_PREC complex *)in + j) = 
00337         ((FFTW_PREC)(*(real_data + j)) + 
00338          (FFTW_PREC)(*(imag_data + j)) * I);
00339     }
00340     }
00341 
00342     /* Now we proceed with the specific processing of each sub-case */
00343 
00344     if (imag_in == NULL) {
00345     /* Sub-case 1 */
00346 
00347     rfftwnd_plan rp = 
00348         rfftw2d_create_plan(nx, ny, FFTW_FORWARD, FFTW_ESTIMATE);
00349         
00350 
00351     /*
00352      * If the input pixel type is coherent with the precision
00353      * of the FFTW library available, the pixel buffer can be 
00354      * directly handled by FFTW. Otherwise, it must be casted.
00355      */
00356 
00357 #if FFTW_PREC_CODE == IN_TYPE_CODE
00358     const IN_TYPE * real_data = 
00359         IN_TYPE_CONST_ADD(cpl_image_get_data)(real_in);
00360         
00361     rfftwnd_real_to_complex(rp, 1, 
00362                 (fftw_real *)real_data, 1, nx * ny,
00363                 out,                    1, nx * ny);
00364 #else
00365     const cpl_image * fprec_img = cpl_image_cast(real_in, CPL_TYPE_ADDED);
00366 
00367     const FFTW_PREC * real_data = 
00368         (FFTW_PREC *)cpl_image_get_data_const(fprec_img);
00369         
00370     rfftwnd_real_to_complex(rp, 1, 
00371                 (fftw_real *)real_data, 1, nx * ny,
00372                 out,                    1, nx * ny);
00373     cpl_image_delete((cpl_image *)fprec_img);
00374 #endif
00375 
00376     rfftwnd_destroy_plan(rp);    
00377     } else if (imag_out == NULL) {
00378     /* Sub-case 2 */
00379         
00380     rfftwnd_plan rp = 
00381         rfftw2d_create_plan(nx, ny, FFTW_BACKWARD, FFTW_ESTIMATE);
00382         
00383     OUT_TYPE * real_outdata = 
00384         OUT_TYPE_ADD(cpl_image_get_data)(real_out);
00385 
00386     /*
00387      * Similarly to sub-case 2, the output pixel type can only
00388      * by directly written by FFTW if it is coherent with the
00389      * precision of the FFTW library available. Otherwise a new
00390      * pixel buffer of the appropriate precision is created and
00391      * copied to the real output buffer after computation.
00392      */
00393 
00394 #if FFTW_PREC_CODE == OUT_TYPE_CODE
00395     rfftwnd_complex_to_real(rp, 1,
00396                 in,                        1, nx * ny,
00397                 (fftw_real *)real_outdata, 1, nx * ny);
00398 #else
00399     cpl_image * fprec_img = 
00400         cpl_image_new(nx, ny, CPL_TYPE_ADDED);
00401     FFTW_PREC * fp_outdata = 
00402         (FFTW_PREC *)cpl_image_get_data(fprec_img);
00403                   
00404     rfftwnd_complex_to_real(rp, 1, 
00405                 in,                      1, nx * ny,
00406                 (fftw_real *)fp_outdata, 1, nx * ny);
00407 
00408     for(j = 0; j < nx * ny; j++) {
00409         *(real_outdata + j) = (OUT_TYPE)(*(fp_outdata + j));
00410     }
00411 
00412     cpl_image_delete(fprec_img);
00413 #endif
00414 
00415     rfftwnd_destroy_plan(rp);
00416     } else {
00417     /* Sub-case 3 */
00418 
00419     fftwnd_plan p =
00420         fftw2d_create_plan_specific(nx, ny, fftwmode, FFTW_ESTIMATE |
00421                     FFTW_IN_PLACE,
00422                     in, 1, out, 1);
00423     fftwnd(p, 1, in, 1, nx * ny, out, 1, nx * ny);
00424     fftwnd_destroy_plan(p);
00425     }
00426 
00427     /* Output preparation for sub-cases 1 and 3 */
00428     if (imag_out != NULL && imag_in != NULL) {
00429         
00430     OUT_TYPE * real_outdata = 
00431         OUT_TYPE_ADD(cpl_image_get_data)(real_out);
00432 
00433     OUT_TYPE * imag_outdata =
00434         OUT_TYPE_ADD(cpl_image_get_data)(imag_out);
00435 
00436     for(j = 0; j < nx * ny; j++) {
00437         *(real_outdata + j) = 
00438         (OUT_TYPE)((out + j) -> re);
00439 
00440         *(imag_outdata + j) = 
00441         (OUT_TYPE)((out + j) -> im);
00442     }
00443     } else if (imag_in == NULL) {
00444         
00445     OUT_TYPE * real_outdata = 
00446         OUT_TYPE_ADD(cpl_image_get_data)(real_out);
00447 
00448     OUT_TYPE * imag_outdata =
00449         OUT_TYPE_ADD(cpl_image_get_data)(imag_out);
00450 
00451     int k;
00452 
00453     for(j = 0; j < nx; j++) {
00454         *(real_outdata + j) = 
00455         (OUT_TYPE)((out + j) -> re);
00456 
00457         *(imag_outdata + j) = 
00458         (OUT_TYPE)((out + j) -> im);
00459     }
00460 
00461     for(j = 1; j < (ny / 2 + 1); j++) {
00462         for(k = 0; k < nx; k++) {
00463         *(real_outdata + j * nx + k) = 
00464             (OUT_TYPE)((out + j * nx + k) -> re);
00465 
00466         *(real_outdata + (ny - j) * nx + k) = 
00467             (OUT_TYPE)((out + j * nx + k) -> re);
00468 
00469         *(imag_outdata + j * nx + k) = 
00470             (OUT_TYPE)((out + j * nx + k) -> im);
00471 
00472         *(imag_outdata + (ny - j) * nx + k) = 
00473             (OUT_TYPE)((out + j * nx + k) -> im);
00474         }
00475     }
00476     }
00477 
00478     fftw_free(in);
00479     
00480     return cpl_error_get_code();
00481 }
00482 
00483 
00484 

Generated on Fri Jul 3 11:23:57 2009 for NACO Pipeline Reference Manual by  doxygen 1.5.8