MIDAPACK - MIcrowave Data Analysis PACKage 1.0beta
Parallel software tools for high performance CMB DA analysis
toeplitz.h
00001 /*
00002 ** file: toeplitz.h, version 1.0, May 2012  
00003 **
00004 ** Project:  Midapack library, ANR MIDAS'09
00005 ** Purpose:  Provide Toeplitz algebra tools suitable for Cosmic Microwave Background (CMB)
00006 **           data analysis.
00007 **
00008 ** Authors:  Pierre Cargemel, Frederic Dauvergne, Maude Le Jeune, Antoine Rogier, 
00009 **           Radek Stompor (APC, Paris)
00010 **
00011 ** 
00012 ** Copyright (c) 2010-2012 APC CNRS Université Paris Diderot
00013 ** 
00014 ** This program is free software; you can redistribute it and/or modify
00015 ** it under the terms of the GNU General Public License as published by
00016 ** the Free Software Foundation; either version 3 of the License, or
00017 ** (at your option) any later version.
00018 ** This program is distributed in the hope that it will be useful,
00019 ** but WITHOUT ANY WARRANTY; without even the implied warranty of
00020 ** MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
00021 ** GNU General Public License for more details.
00022 ** 
00023 ** You should have received a copy of the GNU General Public License
00024 ** along with this program; if not, see http://www.gnu.org/licenses/gpl.html
00025 **
00026 ** For more information about ANR MIDAS'09 project see
00027 ** http://www.apc.univ-paris7.fr/APC_CS/Recherche/Adamis/MIDAS09/index.html
00028 ** 
00029 ** ACKNOWLEDGMENT: This work has been supported in part by the French National 
00030 ** Research Agency (ANR) through COSINUS program (project MIDAS no. ANR-09-COSI-009).
00031 **
00032 **
00033 */
00034 
00035 #ifndef         TOEPLITZ_H_
00036 #define         TOEPLITZ_H_
00037 
00038 #include <mpi.h>
00039 #include <fftw3.h>
00040 #include <omp.h>
00041 #include <stdlib.h>
00042 #include <stdio.h>
00043 #include <math.h>
00044 #include <string.h>
00045 
00046 //=========================================================================
00047 //Basic functions definition
00048 #define max(a,b) (((a) > (b)) ? (a) : (b))
00049 #define min(a,b) (((a) < (b)) ? (a) : (b))
00050  
00051 //=========================================================================
00052 //Fixed parameters
00053 
00055 
00057 #ifndef VERBOSE
00058 #define VERBOSE 0 
00059 #endif
00060 
00062 
00064 #ifndef MPI_USER_TAG
00065 #define MPI_USER_TAG 123
00066 #endif
00067 
00068 
00069 //=========================================================================
00070 // User routines definition (API)
00071 
00072 //Sequential routines (group 11)
00073 int tpltz_init(int n, int lambda, int *nfft, int *blocksize, fftw_complex **T_fft, double *T, fftw_complex **V_fft, double **V_rfft, fftw_plan *plan_f, fftw_plan *plan_b);
00074 
00075 int tpltz_cleanup(fftw_complex **T_fft, fftw_complex **V_fft, double **V_rfft,fftw_plan *plan_f, fftw_plan *plan_b);
00076 
00077 int stmm_core(double **V, int n, int m, fftw_complex *T_fft, int blocksize, int lambda, fftw_complex *V_fft, double *V_rfft, int nfft, fftw_plan plan_f, fftw_plan plan_b, int flag_offset);
00078 
00079 int stmm(double **V, int n, int m, int id0, int l, fftw_complex *T_fft, int lambda, fftw_complex *V_fft, double *V_rfft, fftw_plan plan_f, fftw_plan plan_b, int blocksize, int nfft);
00080 
00081 int stbmm(double **V, int *n, int m, int nrow, double *T, int nb_blocks_local, int nb_blocks_all, int *lambda, int *idv, int idp, int local_V_size);
00082 
00083 int gstbmm(double **V, int *n, int m, int nrow, double *T, int nb_blocks_local, int nb_blocks_all, int *lambda, int *idv, int id0p, int local_V_size, int *id0gap, int *lgap, int ngap);
00084 
00085 int reset_gaps(double **V, int id0,int local_V_size, int m, int nrow, int *id0gap, int *lgap, int ngap);
00086 
00087 
00088 //Mpi routines (group 12)
00089 #ifdef MPI
00090 int mpi_stmm(double **V, int n, int m, int id0, int l, double *T, int lambda, MPI_Comm comm);
00091 
00092 int mpi_stbmm(double **V, int *n, int m, int nrow, double *T, int nb_blocks_local, int nb_blocks_all, int *lambda, int *idv, int idp, int local_V_size, MPI_Comm comm);
00093 
00094 int mpi_gstbmm(double **V, int *n, int m, int nrow, double *T, int nb_blocks_local, int nb_blocks_all, int *lambda, int *idv, int id0p, int local_V_size, int *id0gap, int *lgap, int ngap, MPI_Comm comm);
00095 
00096 #endif
00097 
00098 
00099 //=========================================================================
00100 // User routines definition
00101 
00102 //Low level routines (group 21)
00103 int optimal_blocksize(int n, int lambda, int bs_flag);
00104 
00105 int fftw_init_omp_threads();
00106 
00107 int rhs_init_fftw(int *nfft, int fft_size, fftw_complex **V_fft, double **V_rfft, fftw_plan *plan_f, fftw_plan *plan_b, int fftw_flag);
00108 
00109 int circ_init_fftw(double *T, int fft_size, int lambda, fftw_complex **T_fft);
00110 
00111 int scmm_direct(int fft_size, fftw_complex *C_fft, int ncol, double *V_rfft, double **CV, fftw_complex *V_fft, fftw_plan plan_f_V, fftw_plan plan_b_CV);
00112 
00113 int scmm_basic(double **V, int blocksize, int m, fftw_complex *C_fft, int lambda, double **CV, fftw_complex *V_fft, double *V_rfft, int nfft, fftw_plan plan_f_V, fftw_plan plan_b_CV);
00114 
00115 int stmm_reshape(double **V, int n, int m, int id0, int l, fftw_complex *T_fft, int lambda, fftw_complex *V_fft, double *V_rfft, fftw_plan plan_f, fftw_plan plan_b, int blocksize, int nfft);
00116 
00117 int build_gappy_blocks(int *n, int m, int nrow, double *T, int nb_blocks_local, int nb_blocks_all, int *lambda, int *idv, int *id0gap, int *lgap, int ngap, int *nb_blocks_gappy_final, double *Tgappy, int *idvgappy, int *ngappy, int *lambdagappy, int flag_param_distmin_fixed);
00118 
00119 
00120 //Internal routines (group 22)
00121 int print_error_message(int error_number, char const *file, int line);
00122 
00123 int copy_block(int ninrow, int nincol, double *Vin, int noutrow, int noutcol, double *Vout, int inrow, int incol, int nblockrow, int nblockcol, int outrow, int outcol, double norm, int set_zero_flag);
00124 
00125 int vect2nfftblock(double *V1, int v1_size, double *V2, int fft_size, int nfft, int lambda);
00126 
00127 int nfftblock2vect(double *V2, int fft_size, int nfft, int lambda, double *V1, int v1_size);
00128 
00129 int get_overlapping_blocks_params(int nbloc, int *idv, int *n, int local_V_size, int nrow, int idp, int *idpnew, int *local_V_size_new, int *nnew, int *ifirstBlock, int *ilastBlock);
00130 
00131 
00132 //=========================================================================
00133 #endif      /* !TOEPLITZ_H_ */
00134 
00135