citcossinci.C

00001 /*
00002  *   Copyright (c) 1999-2001 Eric Gourgoulhon
00003  *
00004  *   This file is part of LORENE.
00005  *
00006  *   LORENE 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  *   LORENE 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 LORENE; if not, write to the Free Software
00018  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00019  *
00020  */
00021 
00022 
00023 char citcossinci_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/citcossinci.C,v 1.2 2004/12/27 14:27:28 j_novak Exp $" ;
00024 
00025 /*
00026  * Transformation inverse cos((2*l+1)*theta) ou sin(2*l*theta) (suivant la 
00027  *  parite de l'indice m en phi) sur le deuxieme indice (theta)
00028  *  d'un tableau 3-D representant une fonction antisymetrique par rapport
00029  *  au plan z=0.
00030  *  Utilise la bibliotheque fftw.
00031  *
00032  * Entree:
00033  * -------
00034  *   int* deg   : tableau du nombre effectif de degres de liberte dans chacune 
00035  *        des 3 dimensions: le nombre de points de collocation
00036  *        en theta est  nt = deg[1] et doit etre de la forme
00037  *          nt = 2*p + 1 
00038  *   int* dimc  : tableau du nombre d'elements de cf dans chacune des trois 
00039  *            dimensions.
00040  *        On doit avoir  dimc[1] >= deg[1] = nt. 
00041  *
00042  *   double* cf :   tableau des coefficients c_l de la fonction definis
00043  *            comme suit (a r et phi fixes)
00044  *
00045  *            pour m pair:
00046  *             f(theta) = som_{l=0}^{nt-2} c_l cos( (2 l+1) theta ) . 
00047  *            pour m impair:
00048  *             f(theta) = som_{l=0}^{nt-2} c_l sin( 2 l theta ) . 
00049  *
00050  *            L'espace memoire correspondant a ce
00051  *                        pointeur doit etre dimc[0]*dimc[1]*dimc[2] et doit 
00052  *            avoir ete alloue avant l'appel a la routine.   
00053  *           Le coefficient c_l (0 <= l <= nt-1) doit etre stoke dans
00054  *           le tableau cf comme suit
00055  *                c_l = cf[ dimc[1]*dimc[2] * j + k + dimc[2] * l ]
00056  *           ou j et k sont les indices correspondant a
00057  *           phi et r respectivement.
00058  *
00059  *   int* dimf  : tableau du nombre d'elements de ff dans chacune des trois 
00060  *            dimensions.
00061  *        On doit avoir  dimf[1] >= deg[1] = nt. 
00062  *
00063  * Sortie:
00064  * -------
00065  *   double* ff : tableau des valeurs de la fonction aux nt points de
00066  *                        de collocation
00067  *
00068  *            theta_l =  pi/2 l/(nt-1)       0 <= l <= nt-1 
00069  *
00070  *            L'espace memoire correspondant a ce
00071  *                        pointeur doit etre dimf[0]*dimf[1]*dimf[2] et doit 
00072  *            avoir ete alloue avant l'appel a la routine.   
00073  *            Les valeurs de la fonction sont stokees
00074  *            dans le tableau ff comme suit
00075  *          f( theta_l ) = ff[ dimf[1]*dimf[2] * j + k + dimf[2] * l ]
00076  *           ou j et k sont les indices correspondant a
00077  *           phi et r respectivement.
00078  *
00079  * NB: Si le pointeur cf est egal a ff, la routine ne travaille que sur un 
00080  *     seul tableau, qui constitue une entree/sortie.
00081  *
00082  */
00083 
00084 /*
00085  * $Id: citcossinci.C,v 1.2 2004/12/27 14:27:28 j_novak Exp $
00086  * $Log: citcossinci.C,v $
00087  * Revision 1.2  2004/12/27 14:27:28  j_novak
00088  * Added forgotten "delete [] t1"
00089  *
00090  * Revision 1.1  2004/12/21 17:06:03  j_novak
00091  * Added all files for using fftw3.
00092  *
00093  * Revision 1.4  2003/01/31 10:31:23  e_gourgoulhon
00094  * Suppressed the directive #include <malloc.h> for malloc is defined
00095  * in <stdlib.h>
00096  *
00097  * Revision 1.3  2002/10/16 14:36:53  j_novak
00098  * Reorganization of #include instructions of standard C++, in order to
00099  * use experimental version 3 of gcc.
00100  *
00101  * Revision 1.2  2002/09/09 13:00:40  e_gourgoulhon
00102  * Modification of declaration of Fortran 77 prototypes for
00103  * a better portability (in particular on IBM AIX systems):
00104  * All Fortran subroutine names are now written F77_* and are
00105  * defined in the new file C++/Include/proto_f77.h.
00106  *
00107  * Revision 1.1.1.1  2001/11/20 15:19:28  e_gourgoulhon
00108  * LORENE
00109  *
00110  * Revision 2.0  1999/02/22  15:42:35  hyc
00111  * *** empty log message ***
00112  *
00113  *
00114  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/FFTW3/citcossinci.C,v 1.2 2004/12/27 14:27:28 j_novak Exp $
00115  *
00116  */
00117 
00118 // headers du C
00119 #include <stdlib.h>
00120 #include <fftw3.h>
00121 
00122 //Lorene prototypes
00123 #include "tbl.h"
00124 
00125 // Prototypage des sous-routines utilisees:
00126 fftw_plan back_fft(int, Tbl*&) ;
00127 double* cheb_ini(const int) ;
00128 double* chebimp_ini(const int ) ;
00129 //*****************************************************************************
00130 
00131 void citcossinci(const int* deg, const int* dimc, double* cf, const int* dimf,
00132            double* ff)
00133 {
00134 
00135 int i, j, k ;
00136 
00137 // Dimensions des tableaux ff et cf  :
00138     int n1f = dimf[0] ;
00139     int n2f = dimf[1] ;
00140     int n3f = dimf[2] ;
00141     int n1c = dimc[0] ;
00142     int n2c = dimc[1] ;
00143     int n3c = dimc[2] ;
00144 
00145 // Nombres de degres de liberte en theta :    
00146     int nt = deg[1] ;
00147     
00148 // Tests de dimension:
00149     if (nt > n2f) {
00150     cout << "citcossinci: nt > n2f : nt = " << nt << " ,  n2f = " 
00151     << n2f << endl ;
00152     abort () ;
00153     }
00154     if (nt > n2c) {
00155     cout << "citcossinci: nt > n2c : nt = " << nt << " ,  n2c = " 
00156     << n2c << endl ;
00157     abort () ;
00158     }
00159     if (n1c > n1f) {
00160     cout << "citcossinci: n1c > n1f : n1c = " << n1c << " ,  n1f = " 
00161     << n1f << endl ;
00162     abort () ;
00163     }
00164     if (n3c > n3f) {
00165     cout << "citcossinci: n3c > n3f : n3c = " << n3c << " ,  n3f = " 
00166     << n3f << endl ;
00167     abort () ;
00168     }
00169 
00170 // Nombre de points pour la FFT:
00171     int nm1 = nt - 1;
00172     int nm1s2 = nm1 / 2;
00173 
00174 // Recherche des tables pour la FFT:
00175     Tbl* pg = 0x0 ;
00176     fftw_plan p = back_fft(nm1, pg) ;
00177     Tbl& g = *pg ;
00178     double* t1 = new double[nt] ;
00179 
00180 // Recherche de la table des sin(psi) :
00181     double* sinp = cheb_ini(nt);    
00182     
00183 // Recherche de la table des points de collocations x_k = cos(theta_{nt-1-k}):
00184     double* x = chebimp_ini(nt) ;   
00185 
00186 // boucle sur phi et r (les boucles vont resp. de 0 a dimf[0]-1
00187 //           et 0 a dimf[2])
00188 
00189     int n2n3f = n2f * n3f ;
00190     int n2n3c = n2c * n3c ;
00191 
00192 //=======================================================================
00193 //              Cas m pair
00194 //=======================================================================
00195 
00196     j = 0 ;
00197     
00198     while (j<n1f-1) {   //le dernier coef en phi n'est pas considere
00199             // (car nul)
00200 
00201 //-----------------------------------------------------------------------
00202 //  partie cos(m phi) avec m pair : transformation cos((2 l+1) theta) inverse
00203 //-----------------------------------------------------------------------
00204 
00205     for (k=0; k<n3c; k++) {
00206 
00207         int i0 = n2n3c * j + k ; // indice de depart 
00208         double* cf0 = cf + i0 ;    // tableau des donnees a transformer
00209 
00210         i0 = n2n3f * j + k ; // indice de depart 
00211         double* ff0 = ff + i0 ;    // tableau resultat
00212          
00213 // Calcul des coefficients du developpement en T_{2i}(x) de la fonction
00214 //  h(x) := x f(x)  (x=cos(theta)) a partir des coefficients de f 
00215 // (resultat stoke dans le tableau t1 :
00216     t1[0] = .5 * cf0[0] ; 
00217     for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[n3c*i] + cf0[n3c*(i-1)] ) ;
00218     t1[nm1] = .5 * cf0[n3c*(nt-2)] ;
00219 
00220 /*
00221  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
00222  *     reliee a theta par psi = 2 theta   et F(psi) = f(theta(psi)).  
00223  */
00224  
00225 // Calcul des coefficients de Fourier de la fonction 
00226 //   G(psi) = F+(psi) + F_(psi) sin(psi)
00227 // en fonction des coefficients en cos(2l theta) de f:
00228 
00229 // Coefficients impairs de G
00230 //--------------------------
00231  
00232         double c1 = t1[1] ;
00233 
00234             double som = 0;
00235         ff0[n3f] = 0 ;
00236             for ( i = 3; i < nt; i += 2 ) {
00237             ff0[ n3f*i ] = t1[i] - c1 ;
00238         som += ff0[ n3f*i ] ;
00239             }   
00240 
00241 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
00242         double fmoins0 = nm1s2 * c1 + som ;
00243 
00244 // Coef. impairs de G
00245 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
00246 //     donnait exactement les coef. des sinus, ce facteur serait -0.5.
00247             for ( i = 3; i < nt; i += 2 ) {
00248         g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
00249             }
00250 
00251 
00252 // Coefficients pairs de G
00253 //------------------------
00254 //  Ces coefficients sont egaux aux coefficients pairs du developpement de
00255 //   f.
00256 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
00257 //     donnait exactement les coef. des cosinus, ce facteur serait 1.
00258 
00259         g.set(0) = t1[0] ;
00260             for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * t1[2*i] ; 
00261             g.set(nm1s2) = t1[nm1] ;
00262 
00263 // Transformation de Fourier inverse de G 
00264 //---------------------------------------
00265 
00266 // FFT inverse
00267         fftw_execute(p) ;
00268 
00269 // Valeurs de f deduites de celles de G
00270 //-------------------------------------
00271 
00272             for ( i = 1; i < nm1s2 ; i++ ) {
00273 // ... indice du pt symetrique de psi par rapport a pi/2:
00274         int isym = nm1 - i ; 
00275     
00276         double fp = 0.5 * ( g(i) + g(isym) ) ;
00277         double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
00278         ff0[ n3f*i ] =  ( fp + fm ) / x[isym] ;
00279         ff0[ n3f*isym ] = ( fp - fm ) / x[i] ;
00280             }
00281     
00282 //... cas particuliers:
00283         ff0[0] = g(0) + fmoins0 ;
00284         ff0[ n3f*nm1 ] = 0 ;
00285         ff0[ n3f*nm1s2 ] = g(nm1s2) / x[nm1s2] ;
00286 
00287     }   // fin de la boucle sur r 
00288 
00289 //-----------------------------------------------------------------------
00290 //  partie sin(m phi) avec m pair : transformation cos((2 l+1) theta) inverse
00291 //-----------------------------------------------------------------------
00292 
00293     j++ ;
00294 
00295     if ( (j != 1) && (j != n1f-1 ) ) {  
00296 //  on effectue le calcul seulement dans les cas ou les coef en phi ne sont 
00297 //  pas nuls 
00298 
00299     for (k=0; k<n3c; k++) {
00300 
00301         int i0 = n2n3c * j + k ; // indice de depart 
00302         double* cf0 = cf + i0 ;    // tableau des donnees a transformer
00303 
00304         i0 = n2n3f * j + k ; // indice de depart 
00305         double* ff0 = ff + i0 ;    // tableau resultat
00306          
00307 // Calcul des coefficients du developpement en T_{2i}(x) de la fonction
00308 //  h(x) := x f(x)  (x=cos(theta)) a partir des coefficients de f 
00309 // (resultat stoke dans le tableau t1 :
00310     t1[0] = .5 * cf0[0] ;
00311     for (i=1; i<nm1; i++) t1[i] = .5 * ( cf0[n3c*i] + cf0[n3c*(i-1)] ) ;
00312     t1[nm1] = .5 * cf0[n3c*(nt-2)] ;
00313 
00314 /*
00315  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
00316  *     reliee a theta par psi = 2 theta   et F(psi) = f(theta(psi)).  
00317  */
00318  
00319 // Calcul des coefficients de Fourier de la fonction 
00320 //   G(psi) = F+(psi) + F_(psi) sin(psi)
00321 // en fonction des coefficients en cos(2l theta) de f:
00322 
00323 // Coefficients impairs de G
00324 //--------------------------
00325  
00326         double c1 = t1[1] ;
00327 
00328             double som = 0;
00329         ff0[n3f] = 0 ;
00330             for ( i = 3; i < nt; i += 2 ) {
00331             ff0[ n3f*i ] = t1[i] - c1 ;
00332         som += ff0[ n3f*i ] ;
00333             }   
00334 
00335 // Valeur en psi=0 de la partie antisymetrique de F, F_ :
00336         double fmoins0 = nm1s2 * c1 + som ;
00337 
00338 // Coef. impairs de G
00339 // NB: le facteur 0.25 est du a la normalisation de fftw; si fftw
00340 //     donnait exactement les coef. des sinus, ce facteur serait -0.5.
00341             for ( i = 3; i < nt; i += 2 ) {
00342         g.set(nm1-i/2) = 0.25 * ( ff0[ n3f*i ] - ff0[ n3f*(i-2) ] ) ;
00343             }
00344 
00345 
00346 // Coefficients pairs de G
00347 //------------------------
00348 //  Ces coefficients sont egaux aux coefficients pairs du developpement de
00349 //   f.
00350 // NB: le facteur 0.5 est du a la normalisation de fftw; si fftw
00351 //     donnait exactement les coef. des cosinus, ce facteur serait 1.
00352 
00353         g.set(0) = t1[0] ;
00354             for (i=1; i<nm1s2; i++ ) g.set(i) = 0.5 * t1[2*i] ; 
00355             g.set(nm1s2) = t1[nm1] ;
00356 
00357 // Transformation de Fourier inverse de G 
00358 //---------------------------------------
00359 
00360 // FFT inverse
00361         fftw_execute(p) ;
00362 
00363 // Valeurs de f deduites de celles de G
00364 //-------------------------------------
00365 
00366             for ( i = 1; i < nm1s2 ; i++ ) {
00367 // ... indice du pt symetrique de psi par rapport a pi/2:
00368         int isym = nm1 - i ; 
00369     
00370         double fp = 0.5 * ( g(i) + g(isym) ) ;
00371         double fm = 0.5 * ( g(i) - g(isym) ) / sinp[i] ;
00372         ff0[ n3f*i ] =  ( fp + fm ) / x[isym] ;
00373         ff0[ n3f*isym ] = ( fp - fm ) / x[i] ;
00374             }
00375     
00376 //... cas particuliers:
00377         ff0[0] = g(0) + fmoins0 ;
00378         ff0[ n3f*nm1 ] = 0 ;
00379         ff0[ n3f*nm1s2 ] = g(nm1s2) / x[nm1s2] ;
00380 
00381 
00382     }   // fin de la boucle sur r 
00383 
00384     }    // fin du cas ou le calcul etait necessaire (i.e. ou les
00385          // coef en phi n'etaient pas nuls)
00386 
00387 // On passe au cas m pair suivant:
00388     j+=3 ;
00389 
00390    }    // fin de la boucle sur les cas m pair
00391 
00392 //##
00393     if (n1f<=3) {
00394       delete [] t1 ;
00395       return ;
00396     }
00397     
00398 //=======================================================================
00399 //              Cas m impair
00400 //=======================================================================
00401 
00402     j = 2 ;
00403     
00404     while (j<n1f-1) {   //le dernier coef en phi n'est pas considere
00405             // (car nul)
00406 
00407 //--------------------------------------------------------------------------
00408 //  partie cos(m phi) avec m impair : transformation sin((2 l) theta) inv.
00409 //--------------------------------------------------------------------------
00410 
00411     for (k=0; k<n3c; k++) {
00412 
00413         int i0 = n2n3c * j + k ; // indice de depart 
00414         double* cf0 = cf + i0 ;    // tableau des donnees a transformer
00415 
00416         i0 = n2n3f * j + k ; // indice de depart 
00417         double* ff0 = ff + i0 ;    // tableau resultat
00418          
00419 
00420 /*
00421  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
00422  *     reliee a theta par psi = 2 theta  et F(psi) = f(theta(psi)).  
00423  */
00424  
00425 
00426 // Calcul des coefficients de Fourier de la fonction 
00427 //   G(psi) = F+(psi) + F_(psi) sin(psi)
00428 // en fonction des coefficients en sin(2l theta) de f:
00429 
00430 // Coefficients en sinus de G
00431 //---------------------------
00432 //  Ces coefficients sont egaux aux coefficients pairs du developmt. en
00433 //  sin(2l theta) de f (le facteur -.5 vient de la normalisation
00434 //  de fftw: si fftw donnait reellement les coefficients en sinus,
00435 //  il faudrait le remplacer par un +1) :
00436             for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = - .5 * cf0[n3c*i] ;  
00437 
00438 // Coefficients en cosinus de G
00439 //-----------------------------
00440 // Ces coefficients se deduisent des coefficients impairs du developmt. 
00441 // en sin(2l theta) de f (le facteur +.25 vient de la normalisation
00442 //  de fftw: si fftw donnait reellement les coefficients en cosinus,
00443 //  il faudrait le remplacer par un +.5)
00444 
00445         g.set(0) = .5 * cf0[n3c] ;
00446             for ( i = 1; i < nm1s2; i++ ) {
00447         g.set(i) = .25 * ( cf0[ n3c*(2*i+1) ] - cf0[ n3c*(2*i-1) ] ) ;
00448             }
00449         g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
00450 
00451 
00452 // Transformation de Fourier inverse de G 
00453 //---------------------------------------
00454 
00455 // FFT inverse
00456         fftw_execute(p) ;
00457 
00458 // Valeurs de f deduites de celles de G
00459 //-------------------------------------
00460 
00461             for ( i = 1; i < nm1s2 ; i++ ) {
00462 // ... indice du pt symetrique de psi par rapport a pi/2:
00463         int isym = nm1 - i ; 
00464     
00465         double fp = 0.5 * ( g(i) + g(isym) ) / sinp[i] ;
00466         double fm = 0.5 * ( g(i) - g(isym) ) ;
00467         ff0[ n3f*i ] = fp + fm ;
00468         ff0[ n3f*isym ] = fp - fm ;
00469             }
00470     
00471 //... cas particuliers:
00472         ff0[0] = 0. ;
00473         ff0[ n3f*nm1 ] = -2. * g(0) ;
00474         ff0[ n3f*nm1s2 ] = g(nm1s2) ;
00475 
00476     }   // fin de la boucle sur r 
00477 
00478 
00479 //--------------------------------------------------------------------------
00480 //  partie sin(m phi) avec m impair : transformation sin( (2 l) theta) inv.
00481 //--------------------------------------------------------------------------
00482 
00483     j++ ;
00484 
00485     if ( j != n1f-1  ) {  
00486 //  on effectue le calcul seulement dans les cas ou les coef en phi ne sont 
00487 //  pas nuls 
00488 
00489     for (k=0; k<n3c; k++) {
00490 
00491         int i0 = n2n3c * j + k ; // indice de depart 
00492         double* cf0 = cf + i0 ;    // tableau des donnees a transformer
00493 
00494         i0 = n2n3f * j + k ; // indice de depart 
00495         double* ff0 = ff + i0 ;    // tableau resultat
00496          
00497 
00498 /*
00499  * NB: dans les commentaires qui suivent, psi designe la variable de [0, pi]
00500  *     reliee a theta par psi = 2 theta  et F(psi) = f(theta(psi)).  
00501  */
00502  
00503 
00504 // Calcul des coefficients de Fourier de la fonction 
00505 //   G(psi) = F+(psi) + F_(psi) sin(psi)
00506 // en fonction des coefficients en sin(2l theta) de f:
00507 
00508 // Coefficients en sinus de G
00509 //---------------------------
00510 //  Ces coefficients sont egaux aux coefficients pairs du developmt. en
00511 //  sin(2l theta) de f (le facteur -.5 vient de la normalisation
00512 //  de fftw: si fftw donnait reellement les coefficients en sinus,
00513 //  il faudrait le remplacer par un +1) :
00514             for (i=2; i<nm1; i += 2 ) g.set(nm1-i/2) = - .5 * cf0[n3c*i] ;  
00515 
00516 // Coefficients en cosinus de G
00517 //-----------------------------
00518 // Ces coefficients se deduisent des coefficients impairs du developmt. 
00519 // en sin(2l theta) de f (le facteur +.25 vient de la normalisation
00520 //  de fftw: si fftw donnait reellement les coefficients en cosinus,
00521 //  il faudrait le remplacer par un +.5)
00522 
00523         g.set(0) = .5 * cf0[n3c] ;
00524             for ( i = 1; i < nm1s2; i++ ) {
00525         g.set(i) = .25 * ( cf0[ n3c*(2*i+1) ] - cf0[ n3c*(2*i-1) ] ) ;
00526             }
00527         g.set(nm1s2) = - .5 * cf0[ n3c*(nt-2) ] ;
00528 
00529 
00530 // Transformation de Fourier inverse de G 
00531 //---------------------------------------
00532 
00533 // FFT inverse
00534         fftw_execute(p) ;
00535 
00536 // Valeurs de f deduites de celles de G
00537 //-------------------------------------
00538 
00539             for ( i = 1; i < nm1s2 ; i++ ) {
00540 // ... indice du pt symetrique de psi par rapport a pi/2:
00541         int isym = nm1 - i ; 
00542     
00543         double fp = 0.5 * ( g(i) + g(isym) ) / sinp[i] ;
00544         double fm = 0.5 * ( g(i) - g(isym) ) ;
00545         ff0[ n3f*i ] = fp + fm ;
00546         ff0[ n3f*isym ] = fp - fm ;
00547             }
00548     
00549 //... cas particuliers:
00550         ff0[0] = 0. ;
00551         ff0[ n3f*nm1 ] = -2. * g(0) ;
00552         ff0[ n3f*nm1s2 ] = g(nm1s2) ;
00553 
00554     }   // fin de la boucle sur r 
00555 
00556 
00557     }    // fin du cas ou le calcul etait necessaire (i.e. ou les
00558          // coef en phi n'etaient pas nuls)
00559 
00560 // On passe au cas m impair suivant:
00561     j+=3 ;
00562 
00563    }    // fin de la boucle sur les cas m impair
00564 
00565     delete [] t1 ;
00566 
00567 }

Generated on Tue Feb 7 01:35:15 2012 for LORENE by  doxygen 1.4.6