File: coef_cheb_fft.C

package info (click to toggle)
lorene 0.0.0~cvs20161116%2Bdfsg-1
  • links: PTS, VCS
  • area: main
  • in suites: bookworm, bullseye, buster, stretch
  • size: 26,444 kB
  • ctags: 13,953
  • sloc: cpp: 212,946; fortran: 21,645; makefile: 1,750; sh: 4
file content (198 lines) | stat: -rw-r--r-- 5,722 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
/*
 * Chebyshev transform via a FFT 
 *
 */

/*
 *   Copyright (c) 2005 Eric Gourgoulhon & Jerome Novak
 *
 *   This file is part of LORENE.
 *
 *   LORENE is free software; you can redistribute it and/or modify
 *   it under the terms of the GNU General Public License as published by
 *   the Free Software Foundation; either version 2 of the License, or
 *   (at your option) any later version.
 *
 *   LORENE is distributed in the hope that it will be useful,
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 *   GNU General Public License for more details.
 *
 *   You should have received a copy of the GNU General Public License
 *   along with LORENE; if not, write to the Free Software
 *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
 *
 */

/*
 * $Id: coef_cheb_fft.C,v 1.3 2014/10/06 15:09:47 j_novak Exp $
 * $Log: coef_cheb_fft.C,v $
 * Revision 1.3  2014/10/06 15:09:47  j_novak
 * Modified #include directives to use c++ syntax.
 *
 * Revision 1.2  2005/11/14 14:12:10  e_gourgoulhon
 * Added include <assert.h>
 *
 * Revision 1.1  2005/11/14 01:56:58  e_gourgoulhon
 * First version
 *
 *
 * $Header: /cvsroot/Lorene/School05/Monday/coef_cheb_fft.C,v 1.3 2014/10/06 15:09:47 j_novak Exp $
 *
 */


#include <iostream>

using namespace std ;

#include <cstdlib>
#include <cmath>
#include <cassert>
#include <fftw3.h>

fftw_plan prepare_fft(int n, double*& pg) {

    static const int nmax = 50 ; //Maximal number of FFT sizes 
    static int nworked = 0 ;
    static double* tab_tab[nmax] ;
    static fftw_plan plan_fft[nmax] ;
    static int nb_fft[nmax] ;

  int index = -1 ;
  for (int i=0; ((i<nworked) && (index<0)); i++) 
    if (nb_fft[i] == n) index = i ; //Has the plan already been estimated?

  if (index <0) { //New plan needed
    index = nworked ;
    if (index >= nmax) {
      cout << "prepare_fft: " << endl ;
      cout << "too many plans!" << endl ;
      abort() ;
    }
        
    double* tab = new double[n] ;
    tab_tab[index] = tab ; 

    plan_fft[index] = 
      fftw_plan_r2r_1d(n, tab, tab, FFTW_R2HC, FFTW_MEASURE) ;
      
    nb_fft[index] = n ;
    nworked++ ;
  }
  assert((index>=0)&&(index<nmax)) ;
  pg = tab_tab[index] ;
  return plan_fft[index] ;
}


double* cheb_ini(int n) {

    const int nmax = 30 ; /* Nombre maximun de dimensions differentes */
    static	double*	table_sin[nmax] ;	/* Tableau des pointeurs sur les tableaux */
    static	int	nwork = 0 ;		/* Nombre de tableaux deja initialises */
    static	int	tbn[nmax] ;		/* Tableau des points deja initialises */

    // Ce nombre de points a-t-il deja ete utilise ?
    int indice = -1 ;
    for (int i=0 ; i < nwork ; i++ ) {
	    if ( tbn[i] == n ) indice = i ;
	}

    // Initialisation
    if (indice == -1) {		    /* Il faut une nouvelle initialisation */
	    if ( nwork >= nmax ) {
	        cout << "cheb_ini : nwork >= nmax !" << endl ; 
	        abort() ; 
	    }
	    indice = nwork ; nwork++ ; tbn[indice] = n ;

	    int nm1s2 = (n-1) / 2 ;  		
	    table_sin[indice] = new double[nm1s2] ; 

	    double xx = M_PI / double(n-1);
	    for (int i = 0; i < nm1s2 ; i++ ) {
	        table_sin[indice][i] = sin( xx * i );
	    }
	}

    return table_sin[indice] ;
}


void coef_cheb_fft(int nr, const double* ff, double* cf) {

    // Nombre de points pour la FFT:
    int nm1 = nr - 1;
    int nm1s2 = nm1 / 2;

    // Recherche des tables pour la FFT:
    double* pg = 0x0 ;
    fftw_plan p = prepare_fft(nm1, pg) ;
    double* g = pg ;

    // Recherche de la table des sin(psi) :
    double* sinp = cheb_ini(nr);		

/*
 * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
 *     reliee a x par  x = - cos(psi)   et F(psi) = f(x(psi)).  
 */
 
    // Valeur en psi=0 de la partie antisymetrique de F, F_ :
    double fmoins0 = 0.5 * ( ff[0] - ff[nm1] );

// Fonction G(psi) = F+(psi) + F_(psi) sin(psi) 
//---------------------------------------------
    	    for (int i = 1; i < nm1s2 ; i++ ) {
// ... indice du pt symetrique de psi par rapport a pi/2:
    int isym = nm1 - i ; 
// ... F+(psi)	
    double fp = 0.5 * ( ff[i] + ff[isym] ) ;	
// ... F_(psi) sin(psi)
    double fms = 0.5 * ( ff[i] - ff[isym] ) * sinp[i] ; 
    g[i] = fp + fms ;
    g[isym] = fp - fms ;
    	    }
//... cas particuliers:
    g[0] = 0.5 * ( ff[0] + ff[nm1] );
    g[nm1s2] = ff[nm1s2];

// Developpement de G en series de Fourier par une FFT
//----------------------------------------------------

    fftw_execute(p) ;

// Coefficients pairs du developmt. de Tchebyshev de f
//----------------------------------------------------
//  Ces coefficients sont egaux aux coefficients en cosinus du developpement
//  de G en series de Fourier (le facteur 2/nm1 vient de la normalisation
//  de fftw) :

	    double fac = 2./double(nm1) ;
	    cf[0] = g[0] / double(nm1) ;
	    for (int i=2; i<nm1; i += 2) cf[i] = fac*g[i/2] ;
	    cf[nm1] = g[nm1s2] / double(nm1) ;

// Coefficients impairs du developmt. de Tchebyshev de f
//------------------------------------------------------
// 1. Coef. c'_k (recurrence amorcee a partir de zero):
//    NB: Le 4/nm1 en facteur de g(i) est du a la normalisation de fftw
//  (si fftw donnait reellement les coef. en sinus, il faudrait le
//   remplacer par un +2.) 
	    fac *= -2. ;
    	    cf[1] = 0 ;
    	    double som = 0;
    	    for (int i = 3; i < nr; i += 2 ) {
	      cf[i] = cf[i-2] + fac * g[nm1-i/2] ;
	    	som += cf[i] ;
    	    }

// 2. Calcul de c_1 :
	    double c1 = - ( fmoins0 + som ) / nm1s2 ;

// 3. Coef. c_k avec k impair:	
    	    cf[1] = c1 ;
    	    for (int i = 3; i < nr; i += 2 ) cf[i] += c1 ;

}