// @(#)root/fft:$Id$
// Author: Anna Kreshuk   07/4/2006

/*************************************************************************
 * Copyright (C) 1995-2006, Rene Brun and Fons Rademakers.               *
 * All rights reserved.                                                  *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/LICENSE.                         *
 * For the list of contributors see $ROOTSYS/README/CREDITS.             *
 *************************************************************************/

//////////////////////////////////////////////////////////////////////////
//
// TFFTReal
// One of the interface classes to the FFTW package, can be used directly
// or via the TVirtualFFT class. Only the basic interface of FFTW is implemented.
//
// Computes transforms called r2r in FFTW manual:
// - transforms of real input and output in "halfcomplex" format i.e.
//   real and imaginary parts for a transform of size n stored as
//   (r0, r1, r2, ..., rn/2, i(n+1)/2-1, ..., i2, i1)
// - discrete Hartley transform
// - sine and cosine transforms (DCT-I,II,III,IV and DST-I,II,III,IV)
// For the detailed information on the computed
// transforms please refer to the FFTW manual, chapter "What FFTW really computes".
//
// How to use it:
// 1) Create an instance of TFFTReal - this will allocate input and output
//    arrays (unless an in-place transform is specified)
// 2) Run the Init() function with the desired flags and settings (see function
//    comments for possible kind parameters)
// 3) Set the data (via SetPoints()or SetPoint() functions)
// 4) Run the Transform() function
// 5) Get the output (via GetPoints() or GetPoint() functions)
// 6) Repeat steps 3)-5) as needed
// For a transform of the same size, but of different kind (or with different flags),
// rerun the Init() function and continue with steps 3)-5)
//
// NOTE: 1) running Init() function will overwrite the input array! Don't set any data
//          before running the Init() function!
//       2) FFTW computes unnormalized transform, so doing a transform followed by
//          its inverse will lead to the original array scaled BY:
//          - transform size (N) for R2HC, HC2R, DHT transforms
//          - 2*(N-1) for DCT-I (REDFT00)
//          - 2*(N+1) for DST-I (RODFT00)
//          - 2*N for the remaining transforms
// Transform inverses:
// R2HC<-->HC2R
// DHT<-->DHT
// DCT-I<-->DCT-I
// DCT-II<-->DCT-III
// DCT-IV<-->DCT-IV
// DST-I<-->DST-I
// DST-II<-->DST-III
// DST-IV<-->DST-IV
//
//////////////////////////////////////////////////////////////////////////

#include "TFFTReal.h"
#include "fftw3.h"

ClassImp(TFFTReal)

//_____________________________________________________________________________
TFFTReal::TFFTReal()
{
//default

   fIn    = 0;
   fOut   = 0;
   fPlan  = 0;
   fN     = 0;
   fKind  = 0;
   fFlags = 0;
   fNdim = 0;
   fTotalSize = 0;
}

//_____________________________________________________________________________
TFFTReal::TFFTReal(Int_t n, Bool_t inPlace)
{
//For 1d transforms
//n here is the physical size of the transform (see FFTW manual for more details)

   fIn = fftw_malloc(sizeof(Double_t)*n);
   if (inPlace) fOut = 0;
   else fOut = fftw_malloc(sizeof(Double_t)*n);

   fPlan = 0;
   fNdim = 1;
   fN = new Int_t[1];
   fN[0] = n;
   fKind = 0;
   fTotalSize = n;
   fFlags = 0;
}

//_____________________________________________________________________________
TFFTReal::TFFTReal(Int_t ndim, Int_t *n, Bool_t inPlace)
{
//For multidimensional transforms
//1st parameter is the # of dimensions,
//2nd is the sizes (physical) of the transform in each dimension

   fTotalSize = 1;
   fNdim = ndim;
   fN = new Int_t[ndim];
   fKind = 0;
   fPlan = 0;
   fFlags = 0;
   for (Int_t i=0; i<ndim; i++){
      fTotalSize*=n[i];
      fN[i] = n[i];
   }
   fIn = fftw_malloc(sizeof(Double_t)*fTotalSize);
   if (!inPlace)
      fOut = fftw_malloc(sizeof(Double_t)*fTotalSize);
   else
      fOut = 0;
}

//_____________________________________________________________________________
TFFTReal::~TFFTReal()
{
//clean-up

   fftw_destroy_plan((fftw_plan)fPlan);
   fPlan = 0;
   fftw_free(fIn);
   fIn = 0;
   if (fOut){
      fftw_free(fOut);
   }
   fOut = 0;
   if (fN)
      delete [] fN;
   fN = 0;
   if (fKind)
      fftw_free((fftw_r2r_kind*)fKind);
   fKind = 0;
}

//_____________________________________________________________________________
void TFFTReal::Init( Option_t* flags,Int_t /*sign*/, const Int_t *kind)
{
//Creates the fftw-plan
//
//NOTE:  input and output arrays are overwritten during initialisation,
//       so don't set any points, before running this function!!!!!
//
//1st parameter:
//  Possible flag_options:
//  "ES" (from "estimate") - no time in preparing the transform, but probably sub-optimal
//       performance
//  "M" (from "measure") - some time spend in finding the optimal way to do the transform
//  "P" (from "patient") - more time spend in finding the optimal way to do the transform
//  "EX" (from "exhaustive") - the most optimal way is found
//  This option should be chosen depending on how many transforms of the same size and
//  type are going to be done. Planning is only done once, for the first transform of this
//  size and type.
//2nd parameter is dummy and doesn't need to be specified
//3rd parameter- transform kind for each dimension
//     4 different kinds of sine and cosine transforms are available
//     DCT-I   - kind=0
//     DCT-II  - kind=1
//     DCT-III - kind=2
//     DCT-IV  - kind=3
//     DST-I   - kind=4
//     DST-II  - kind=5
//     DSTIII  - kind=6
//     DSTIV   - kind=7

   if (fPlan)
      fftw_destroy_plan((fftw_plan)fPlan);
   fPlan = 0;

   if (!fKind)
      fKind = (fftw_r2r_kind*)fftw_malloc(sizeof(fftw_r2r_kind)*fNdim);

   if (MapOptions(kind)){
      if (fOut)
         fPlan = (void*)fftw_plan_r2r(fNdim, fN, (Double_t*)fIn, (Double_t*)fOut, (fftw_r2r_kind*)fKind, MapFlag(flags));
      else
         fPlan = (void*)fftw_plan_r2r(fNdim, fN, (Double_t*)fIn, (Double_t*)fIn, (fftw_r2r_kind*)fKind, MapFlag(flags));
      fFlags = flags;
   }
}

//_____________________________________________________________________________
void TFFTReal::Transform()
{
//Computes the transform, specified in Init() function

   if (fPlan)
      fftw_execute((fftw_plan)fPlan);
   else {
      Error("Transform", "transform hasn't been initialised");
      return;
   }
}

//_____________________________________________________________________________
Option_t *TFFTReal::GetType() const
{
//Returns the type of the transform

   if (!fKind) {
      Error("GetType", "Type not defined yet (kind not set)");
      return "";
   }
   if (((fftw_r2r_kind*)fKind)[0]==FFTW_R2HC) return "R2HC";
   if (((fftw_r2r_kind*)fKind)[0]==FFTW_HC2R) return "HC2R";
   if (((fftw_r2r_kind*)fKind)[0]==FFTW_DHT) return "DHT";
   else return "R2R";
}

//_____________________________________________________________________________
void TFFTReal::GetPoints(Double_t *data, Bool_t fromInput) const
{
//Copies the output (or input) points into the provided array, that should
//be big enough

   const Double_t * array = GetPointsReal(fromInput);
   if (!array) return;
   std::copy(array, array+fTotalSize, data);
}

//_____________________________________________________________________________
Double_t TFFTReal::GetPointReal(Int_t ipoint, Bool_t fromInput) const
{
//For 1d tranforms. Returns point #ipoint

   if (ipoint<0 || ipoint>fTotalSize){
      Error("GetPointReal", "No such point");
      return 0;
   }
   const Double_t * array = GetPointsReal(fromInput);
   return ( array ) ? array[ipoint] : 0;
}

//_____________________________________________________________________________
Double_t TFFTReal::GetPointReal(const Int_t *ipoint, Bool_t fromInput) const
{
//For multidim.transforms. Returns point #ipoint

   Int_t ireal = ipoint[0];
   for (Int_t i=0; i<fNdim-1; i++)
      ireal=fN[i+1]*ireal + ipoint[i+1];

   const Double_t * array = GetPointsReal(fromInput);
   return ( array ) ? array[ireal] : 0;
}

//_____________________________________________________________________________
void TFFTReal::GetPointComplex(Int_t ipoint, Double_t &re, Double_t &im, Bool_t fromInput) const
{
//Only for input of HC2R and output of R2HC

   const Double_t * array = GetPointsReal(fromInput);
   if (!array) return;
   if ( ( ((fftw_r2r_kind*)fKind)[0]==FFTW_R2HC && !fromInput ) ||
        ( ((fftw_r2r_kind*)fKind)[0]==FFTW_HC2R &&  fromInput ) )
   {
      if (ipoint<fN[0]/2+1){
         re = array[ipoint];
         im = array[fN[0]-ipoint];
      } else {
         re = array[fN[0]-ipoint];
         im = -array[ipoint];
      }
      if ((fN[0]%2)==0 && ipoint==fN[0]/2) im = 0;
   }
}
//_____________________________________________________________________________
void TFFTReal::GetPointComplex(const Int_t *ipoint, Double_t &re, Double_t &im, Bool_t fromInput) const
{
//Only for input of HC2R and output of R2HC and for 1d

   GetPointComplex(ipoint[0], re, im, fromInput);
}

//_____________________________________________________________________________
Double_t* TFFTReal::GetPointsReal(Bool_t fromInput) const
{
//Returns the output (or input) array

   // we have 4 different cases
   // fromInput = false; fOut = !NULL (transformed is not in place) : return fOut
   // fromInput = false; fOut = NULL (transformed is in place) : return fIn
   // fromInput = true; fOut = !NULL :   return fIn
   // fromInput = true; fOut = NULL return an error since input array is overwritten

   if (!fromInput && fOut)
      return (Double_t*)fOut;
   else if (fromInput && !fOut) {
      Error("GetPointsReal","Input array was destroyed");
      return 0;
   }
   //R__ASSERT(fIn);
   return (Double_t*)fIn;
}

//_____________________________________________________________________________
void TFFTReal::SetPoint(Int_t ipoint, Double_t re, Double_t im)
{
   if (ipoint<0 || ipoint>fTotalSize){
      Error("SetPoint", "illegal point index");
      return;
   }
   if (((fftw_r2r_kind*)fKind)[0]==FFTW_HC2R){
      if ((fN[0]%2)==0 && ipoint==fN[0]/2)
         ((Double_t*)fIn)[ipoint] = re;
      else {
         ((Double_t*)fIn)[ipoint] = re;
         ((Double_t*)fIn)[fN[0]-ipoint]=im;
      }
   }
   else
      ((Double_t*)fIn)[ipoint]=re;
}

//_____________________________________________________________________________
void TFFTReal::SetPoint(const Int_t *ipoint, Double_t re, Double_t /*im*/)
{
//Since multidimensional R2HC and HC2R transforms are not supported,
//third parameter is dummy

   Int_t ireal = ipoint[0];
   for (Int_t i=0; i<fNdim-1; i++)
      ireal=fN[i+1]*ireal + ipoint[i+1];
   if (ireal < 0 || ireal >fTotalSize){
      Error("SetPoint", "illegal point index");
      return;
   }
   ((Double_t*)fIn)[ireal]=re;
}

//_____________________________________________________________________________
void TFFTReal::SetPoints(const Double_t *data)
{
//Sets all points

   for (Int_t i=0; i<fTotalSize; i++)
      ((Double_t*)fIn)[i] = data[i];
}

//_____________________________________________________________________________
Int_t TFFTReal::MapOptions(const Int_t *kind)
{
//transfers the r2r_kind parameters to fftw type

   if (kind[0] == 10){
      if (fNdim>1){
         Error("Init", "Multidimensional R2HC transforms are not supported, use R2C interface instead");
         return 0;
      }
      ((fftw_r2r_kind*)fKind)[0] = FFTW_R2HC;
   }
   else if (kind[0] == 11) {
      if (fNdim>1){
         Error("Init", "Multidimensional HC2R transforms are not supported, use C2R interface instead");
         return 0;
      }
      ((fftw_r2r_kind*)fKind)[0] = FFTW_HC2R;
   }
   else if (kind[0] == 12) {
      for (Int_t i=0; i<fNdim; i++)
         ((fftw_r2r_kind*)fKind)[i] = FFTW_DHT;
   }
   else {
      for (Int_t i=0; i<fNdim; i++){
         switch (kind[i]) {
         case 0: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT00;  break;
         case 1: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT01;  break;
         case 2: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT10;  break;
         case 3: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT11;  break;
         case 4: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT00;  break;
         case 5: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT01;  break;
         case 6: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT10;  break;
         case 7: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT11;  break;
         default:
         ((fftw_r2r_kind*)fKind)[i] = FFTW_R2HC; break;
         }
      }
   }
   return 1;
}

//_____________________________________________________________________________
UInt_t TFFTReal::MapFlag(Option_t *flag)
{
//allowed options:
//"ES" - FFTW_ESTIMATE
//"M" - FFTW_MEASURE
//"P" - FFTW_PATIENT
//"EX" - FFTW_EXHAUSTIVE

   TString opt = flag;
   opt.ToUpper();
   if (opt.Contains("ES"))
      return FFTW_ESTIMATE;
   if (opt.Contains("M"))
      return FFTW_MEASURE;
   if (opt.Contains("P"))
      return FFTW_PATIENT;
   if (opt.Contains("EX"))
      return FFTW_EXHAUSTIVE;
   return FFTW_ESTIMATE;
}
 TFFTReal.cxx:1
 TFFTReal.cxx:2
 TFFTReal.cxx:3
 TFFTReal.cxx:4
 TFFTReal.cxx:5
 TFFTReal.cxx:6
 TFFTReal.cxx:7
 TFFTReal.cxx:8
 TFFTReal.cxx:9
 TFFTReal.cxx:10
 TFFTReal.cxx:11
 TFFTReal.cxx:12
 TFFTReal.cxx:13
 TFFTReal.cxx:14
 TFFTReal.cxx:15
 TFFTReal.cxx:16
 TFFTReal.cxx:17
 TFFTReal.cxx:18
 TFFTReal.cxx:19
 TFFTReal.cxx:20
 TFFTReal.cxx:21
 TFFTReal.cxx:22
 TFFTReal.cxx:23
 TFFTReal.cxx:24
 TFFTReal.cxx:25
 TFFTReal.cxx:26
 TFFTReal.cxx:27
 TFFTReal.cxx:28
 TFFTReal.cxx:29
 TFFTReal.cxx:30
 TFFTReal.cxx:31
 TFFTReal.cxx:32
 TFFTReal.cxx:33
 TFFTReal.cxx:34
 TFFTReal.cxx:35
 TFFTReal.cxx:36
 TFFTReal.cxx:37
 TFFTReal.cxx:38
 TFFTReal.cxx:39
 TFFTReal.cxx:40
 TFFTReal.cxx:41
 TFFTReal.cxx:42
 TFFTReal.cxx:43
 TFFTReal.cxx:44
 TFFTReal.cxx:45
 TFFTReal.cxx:46
 TFFTReal.cxx:47
 TFFTReal.cxx:48
 TFFTReal.cxx:49
 TFFTReal.cxx:50
 TFFTReal.cxx:51
 TFFTReal.cxx:52
 TFFTReal.cxx:53
 TFFTReal.cxx:54
 TFFTReal.cxx:55
 TFFTReal.cxx:56
 TFFTReal.cxx:57
 TFFTReal.cxx:58
 TFFTReal.cxx:59
 TFFTReal.cxx:60
 TFFTReal.cxx:61
 TFFTReal.cxx:62
 TFFTReal.cxx:63
 TFFTReal.cxx:64
 TFFTReal.cxx:65
 TFFTReal.cxx:66
 TFFTReal.cxx:67
 TFFTReal.cxx:68
 TFFTReal.cxx:69
 TFFTReal.cxx:70
 TFFTReal.cxx:71
 TFFTReal.cxx:72
 TFFTReal.cxx:73
 TFFTReal.cxx:74
 TFFTReal.cxx:75
 TFFTReal.cxx:76
 TFFTReal.cxx:77
 TFFTReal.cxx:78
 TFFTReal.cxx:79
 TFFTReal.cxx:80
 TFFTReal.cxx:81
 TFFTReal.cxx:82
 TFFTReal.cxx:83
 TFFTReal.cxx:84
 TFFTReal.cxx:85
 TFFTReal.cxx:86
 TFFTReal.cxx:87
 TFFTReal.cxx:88
 TFFTReal.cxx:89
 TFFTReal.cxx:90
 TFFTReal.cxx:91
 TFFTReal.cxx:92
 TFFTReal.cxx:93
 TFFTReal.cxx:94
 TFFTReal.cxx:95
 TFFTReal.cxx:96
 TFFTReal.cxx:97
 TFFTReal.cxx:98
 TFFTReal.cxx:99
 TFFTReal.cxx:100
 TFFTReal.cxx:101
 TFFTReal.cxx:102
 TFFTReal.cxx:103
 TFFTReal.cxx:104
 TFFTReal.cxx:105
 TFFTReal.cxx:106
 TFFTReal.cxx:107
 TFFTReal.cxx:108
 TFFTReal.cxx:109
 TFFTReal.cxx:110
 TFFTReal.cxx:111
 TFFTReal.cxx:112
 TFFTReal.cxx:113
 TFFTReal.cxx:114
 TFFTReal.cxx:115
 TFFTReal.cxx:116
 TFFTReal.cxx:117
 TFFTReal.cxx:118
 TFFTReal.cxx:119
 TFFTReal.cxx:120
 TFFTReal.cxx:121
 TFFTReal.cxx:122
 TFFTReal.cxx:123
 TFFTReal.cxx:124
 TFFTReal.cxx:125
 TFFTReal.cxx:126
 TFFTReal.cxx:127
 TFFTReal.cxx:128
 TFFTReal.cxx:129
 TFFTReal.cxx:130
 TFFTReal.cxx:131
 TFFTReal.cxx:132
 TFFTReal.cxx:133
 TFFTReal.cxx:134
 TFFTReal.cxx:135
 TFFTReal.cxx:136
 TFFTReal.cxx:137
 TFFTReal.cxx:138
 TFFTReal.cxx:139
 TFFTReal.cxx:140
 TFFTReal.cxx:141
 TFFTReal.cxx:142
 TFFTReal.cxx:143
 TFFTReal.cxx:144
 TFFTReal.cxx:145
 TFFTReal.cxx:146
 TFFTReal.cxx:147
 TFFTReal.cxx:148
 TFFTReal.cxx:149
 TFFTReal.cxx:150
 TFFTReal.cxx:151
 TFFTReal.cxx:152
 TFFTReal.cxx:153
 TFFTReal.cxx:154
 TFFTReal.cxx:155
 TFFTReal.cxx:156
 TFFTReal.cxx:157
 TFFTReal.cxx:158
 TFFTReal.cxx:159
 TFFTReal.cxx:160
 TFFTReal.cxx:161
 TFFTReal.cxx:162
 TFFTReal.cxx:163
 TFFTReal.cxx:164
 TFFTReal.cxx:165
 TFFTReal.cxx:166
 TFFTReal.cxx:167
 TFFTReal.cxx:168
 TFFTReal.cxx:169
 TFFTReal.cxx:170
 TFFTReal.cxx:171
 TFFTReal.cxx:172
 TFFTReal.cxx:173
 TFFTReal.cxx:174
 TFFTReal.cxx:175
 TFFTReal.cxx:176
 TFFTReal.cxx:177
 TFFTReal.cxx:178
 TFFTReal.cxx:179
 TFFTReal.cxx:180
 TFFTReal.cxx:181
 TFFTReal.cxx:182
 TFFTReal.cxx:183
 TFFTReal.cxx:184
 TFFTReal.cxx:185
 TFFTReal.cxx:186
 TFFTReal.cxx:187
 TFFTReal.cxx:188
 TFFTReal.cxx:189
 TFFTReal.cxx:190
 TFFTReal.cxx:191
 TFFTReal.cxx:192
 TFFTReal.cxx:193
 TFFTReal.cxx:194
 TFFTReal.cxx:195
 TFFTReal.cxx:196
 TFFTReal.cxx:197
 TFFTReal.cxx:198
 TFFTReal.cxx:199
 TFFTReal.cxx:200
 TFFTReal.cxx:201
 TFFTReal.cxx:202
 TFFTReal.cxx:203
 TFFTReal.cxx:204
 TFFTReal.cxx:205
 TFFTReal.cxx:206
 TFFTReal.cxx:207
 TFFTReal.cxx:208
 TFFTReal.cxx:209
 TFFTReal.cxx:210
 TFFTReal.cxx:211
 TFFTReal.cxx:212
 TFFTReal.cxx:213
 TFFTReal.cxx:214
 TFFTReal.cxx:215
 TFFTReal.cxx:216
 TFFTReal.cxx:217
 TFFTReal.cxx:218
 TFFTReal.cxx:219
 TFFTReal.cxx:220
 TFFTReal.cxx:221
 TFFTReal.cxx:222
 TFFTReal.cxx:223
 TFFTReal.cxx:224
 TFFTReal.cxx:225
 TFFTReal.cxx:226
 TFFTReal.cxx:227
 TFFTReal.cxx:228
 TFFTReal.cxx:229
 TFFTReal.cxx:230
 TFFTReal.cxx:231
 TFFTReal.cxx:232
 TFFTReal.cxx:233
 TFFTReal.cxx:234
 TFFTReal.cxx:235
 TFFTReal.cxx:236
 TFFTReal.cxx:237
 TFFTReal.cxx:238
 TFFTReal.cxx:239
 TFFTReal.cxx:240
 TFFTReal.cxx:241
 TFFTReal.cxx:242
 TFFTReal.cxx:243
 TFFTReal.cxx:244
 TFFTReal.cxx:245
 TFFTReal.cxx:246
 TFFTReal.cxx:247
 TFFTReal.cxx:248
 TFFTReal.cxx:249
 TFFTReal.cxx:250
 TFFTReal.cxx:251
 TFFTReal.cxx:252
 TFFTReal.cxx:253
 TFFTReal.cxx:254
 TFFTReal.cxx:255
 TFFTReal.cxx:256
 TFFTReal.cxx:257
 TFFTReal.cxx:258
 TFFTReal.cxx:259
 TFFTReal.cxx:260
 TFFTReal.cxx:261
 TFFTReal.cxx:262
 TFFTReal.cxx:263
 TFFTReal.cxx:264
 TFFTReal.cxx:265
 TFFTReal.cxx:266
 TFFTReal.cxx:267
 TFFTReal.cxx:268
 TFFTReal.cxx:269
 TFFTReal.cxx:270
 TFFTReal.cxx:271
 TFFTReal.cxx:272
 TFFTReal.cxx:273
 TFFTReal.cxx:274
 TFFTReal.cxx:275
 TFFTReal.cxx:276
 TFFTReal.cxx:277
 TFFTReal.cxx:278
 TFFTReal.cxx:279
 TFFTReal.cxx:280
 TFFTReal.cxx:281
 TFFTReal.cxx:282
 TFFTReal.cxx:283
 TFFTReal.cxx:284
 TFFTReal.cxx:285
 TFFTReal.cxx:286
 TFFTReal.cxx:287
 TFFTReal.cxx:288
 TFFTReal.cxx:289
 TFFTReal.cxx:290
 TFFTReal.cxx:291
 TFFTReal.cxx:292
 TFFTReal.cxx:293
 TFFTReal.cxx:294
 TFFTReal.cxx:295
 TFFTReal.cxx:296
 TFFTReal.cxx:297
 TFFTReal.cxx:298
 TFFTReal.cxx:299
 TFFTReal.cxx:300
 TFFTReal.cxx:301
 TFFTReal.cxx:302
 TFFTReal.cxx:303
 TFFTReal.cxx:304
 TFFTReal.cxx:305
 TFFTReal.cxx:306
 TFFTReal.cxx:307
 TFFTReal.cxx:308
 TFFTReal.cxx:309
 TFFTReal.cxx:310
 TFFTReal.cxx:311
 TFFTReal.cxx:312
 TFFTReal.cxx:313
 TFFTReal.cxx:314
 TFFTReal.cxx:315
 TFFTReal.cxx:316
 TFFTReal.cxx:317
 TFFTReal.cxx:318
 TFFTReal.cxx:319
 TFFTReal.cxx:320
 TFFTReal.cxx:321
 TFFTReal.cxx:322
 TFFTReal.cxx:323
 TFFTReal.cxx:324
 TFFTReal.cxx:325
 TFFTReal.cxx:326
 TFFTReal.cxx:327
 TFFTReal.cxx:328
 TFFTReal.cxx:329
 TFFTReal.cxx:330
 TFFTReal.cxx:331
 TFFTReal.cxx:332
 TFFTReal.cxx:333
 TFFTReal.cxx:334
 TFFTReal.cxx:335
 TFFTReal.cxx:336
 TFFTReal.cxx:337
 TFFTReal.cxx:338
 TFFTReal.cxx:339
 TFFTReal.cxx:340
 TFFTReal.cxx:341
 TFFTReal.cxx:342
 TFFTReal.cxx:343
 TFFTReal.cxx:344
 TFFTReal.cxx:345
 TFFTReal.cxx:346
 TFFTReal.cxx:347
 TFFTReal.cxx:348
 TFFTReal.cxx:349
 TFFTReal.cxx:350
 TFFTReal.cxx:351
 TFFTReal.cxx:352
 TFFTReal.cxx:353
 TFFTReal.cxx:354
 TFFTReal.cxx:355
 TFFTReal.cxx:356
 TFFTReal.cxx:357
 TFFTReal.cxx:358
 TFFTReal.cxx:359
 TFFTReal.cxx:360
 TFFTReal.cxx:361
 TFFTReal.cxx:362
 TFFTReal.cxx:363
 TFFTReal.cxx:364
 TFFTReal.cxx:365
 TFFTReal.cxx:366
 TFFTReal.cxx:367
 TFFTReal.cxx:368
 TFFTReal.cxx:369
 TFFTReal.cxx:370
 TFFTReal.cxx:371
 TFFTReal.cxx:372
 TFFTReal.cxx:373
 TFFTReal.cxx:374
 TFFTReal.cxx:375
 TFFTReal.cxx:376
 TFFTReal.cxx:377
 TFFTReal.cxx:378
 TFFTReal.cxx:379
 TFFTReal.cxx:380
 TFFTReal.cxx:381
 TFFTReal.cxx:382
 TFFTReal.cxx:383
 TFFTReal.cxx:384
 TFFTReal.cxx:385
 TFFTReal.cxx:386
 TFFTReal.cxx:387
 TFFTReal.cxx:388
 TFFTReal.cxx:389
 TFFTReal.cxx:390
 TFFTReal.cxx:391
 TFFTReal.cxx:392
 TFFTReal.cxx:393
 TFFTReal.cxx:394
 TFFTReal.cxx:395
 TFFTReal.cxx:396
 TFFTReal.cxx:397
 TFFTReal.cxx:398
 TFFTReal.cxx:399
 TFFTReal.cxx:400
 TFFTReal.cxx:401
 TFFTReal.cxx:402
 TFFTReal.cxx:403
 TFFTReal.cxx:404
 TFFTReal.cxx:405
 TFFTReal.cxx:406
 TFFTReal.cxx:407
 TFFTReal.cxx:408
 TFFTReal.cxx:409