Public Types | Public Member Functions | Private Attributes | List of all members
Fw2dFFT Class Reference

#include <Fw2dFFT.h>

Public Types

using DataFloat = float
 
using DftFloat = double
 
using Data = Real2dData< DataFloat >
 
using DFT = FftwReal2dDftData< DftFloat >
 
using Index = DFT::Index
 
using IndexArray = DFT::IndexArray
 
using Complex = DFT::Complex
 
using Plan = fftw_plan
 
using PlanMap = std::map< IndexArray, Plan >
 
typedef double FftwComplex[2]
 

Public Member Functions

 Fw2dFFT (Index ndatMax, Index opt)
 
 ~Fw2dFFT ()
 
Index checkDataSize (const IndexArray &nsams) const
 
bool haveForwardPlan (const IndexArray &nsams) const
 
bool haveBackwardPlan (const IndexArray &nsams) const
 
PlanforwardPlan (const IndexArray &nsams)
 
PlanbackwardPlan (const IndexArray &nsams)
 
int fftForward (const Data &dat, DFT &dft, Index logLevel=0)
 
int fftBackward (const DFT &dft, Data &dat, Index logLevel=0)
 
FftwComplexfftwOutData ()
 
DftFloatfloatOutData ()
 

Private Attributes

Index m_ndatMax
 
Index m_flag
 
DftFloatm_inData
 
Complexm_outData
 
PlanMap m_forwardPlans
 
PlanMap m_backwardPlans
 

Detailed Description

Definition at line 30 of file Fw2dFFT.h.

Member Typedef Documentation

Definition at line 40 of file Fw2dFFT.h.

Definition at line 36 of file Fw2dFFT.h.

using Fw2dFFT::DataFloat = float

Definition at line 34 of file Fw2dFFT.h.

Definition at line 37 of file Fw2dFFT.h.

using Fw2dFFT::DftFloat = double

Definition at line 35 of file Fw2dFFT.h.

typedef double Fw2dFFT::FftwComplex[2]

Definition at line 43 of file Fw2dFFT.h.

Definition at line 38 of file Fw2dFFT.h.

Definition at line 39 of file Fw2dFFT.h.

using Fw2dFFT::Plan = fftw_plan

Definition at line 41 of file Fw2dFFT.h.

using Fw2dFFT::PlanMap = std::map<IndexArray, Plan>

Definition at line 42 of file Fw2dFFT.h.

Constructor & Destructor Documentation

Fw2dFFT::Fw2dFFT ( Index  ndatMax,
Index  opt 
)

Definition at line 22 of file Fw2dFFT.cxx.

23 : m_ndatMax(ndatMax),
24  m_flag(opt==2 ? FFTW_PATIENT : opt==1 ? FFTW_MEASURE : FFTW_ESTIMATE),
25  m_inData(reinterpret_cast<DftFloat*>(fftw_malloc(m_ndatMax*sizeof(DftFloat)))),
26  m_outData(reinterpret_cast<Complex*>(fftw_malloc(m_ndatMax*sizeof(DftFloat)))) { }
Index m_flag
Definition: Fw2dFFT.h:91
opt
Definition: train.py:196
DftFloat * m_inData
Definition: Fw2dFFT.h:92
Complex * m_outData
Definition: Fw2dFFT.h:93
double DftFloat
Definition: Fw2dFFT.h:35
Index m_ndatMax
Definition: Fw2dFFT.h:90
Fw2dFFT::~Fw2dFFT ( )

Definition at line 30 of file Fw2dFFT.cxx.

30  {
31  for ( auto& iplan : m_forwardPlans ) fftw_destroy_plan(iplan.second);
32  for ( auto& iplan : m_backwardPlans ) fftw_destroy_plan(iplan.second);
33  fftw_free(m_inData);
34  fftw_free(m_outData);
35 }
PlanMap m_forwardPlans
Definition: Fw2dFFT.h:94
DftFloat * m_inData
Definition: Fw2dFFT.h:92
Complex * m_outData
Definition: Fw2dFFT.h:93
PlanMap m_backwardPlans
Definition: Fw2dFFT.h:95

Member Function Documentation

Fw2dFFT::Plan & Fw2dFFT::backwardPlan ( const IndexArray nsams)

Definition at line 77 of file Fw2dFFT.cxx.

77  {
78  const string myname = "Fw2dFFT::backwardPlan: ";
79  static Plan badplan;
80  if ( checkDataSize(nsams) ) {
81  cout << myname << "Data cannot be accomodated. Maximum data size is " << m_ndatMax << endl;
82  return badplan;
83  }
84  if ( m_backwardPlans.count(nsams) == 0 ) {
85  m_backwardPlans[nsams] = fftw_plan_dft_c2r_2d(nsams[0], nsams[1], fftwOutData(), m_inData, m_flag);
86  }
87  return m_backwardPlans[nsams];
88 }
Index checkDataSize(const IndexArray &nsams) const
Definition: Fw2dFFT.cxx:39
Index m_flag
Definition: Fw2dFFT.h:91
FftwComplex * fftwOutData()
Definition: Fw2dFFT.h:83
DftFloat * m_inData
Definition: Fw2dFFT.h:92
PlanMap m_backwardPlans
Definition: Fw2dFFT.h:95
Index m_ndatMax
Definition: Fw2dFFT.h:90
QTextStream & endl(QTextStream &s)
fftw_plan Plan
Definition: Fw2dFFT.h:41
Fw2dFFT::Index Fw2dFFT::checkDataSize ( const IndexArray nsams) const

Definition at line 39 of file Fw2dFFT.cxx.

39  {
42  if ( ndatIn > m_ndatMax ) return 1;
43  if ( ndatOut > m_ndatMax ) return 2;
44  if ( ndatIn == 0 ) return 3;
45  return 0;
46 }
unsigned int Index
static Index dataSize(const std::array< Index, N > &nsams)
Definition: Real2dDftData.h:34
static Index dftFloatDataSize(const std::array< Index, N > &nsams)
Index m_ndatMax
Definition: Fw2dFFT.h:90
int Fw2dFFT::fftBackward ( const DFT dft,
Data dat,
Index  logLevel = 0 
)

Definition at line 123 of file Fw2dFFT.cxx.

123  {
124  const string myname = "Fw2dFFT::fftBackward: ";
125  if ( ! dft.isValid() ) return 1;
126  auto nsams = dft.nSamples();
127  if ( checkDataSize(nsams) ) {
128  cout << myname << "Sample counts are too large. Maximum data size is " << m_ndatMax << endl;
129  return 2;
130  }
133  if ( dft.normalization().isPower() ) {
134  cout << myname << "ERROR: Power normalization is not (yet) supported." << endl;
135  return 3;
136  }
137  dat.reset(nsams);
138  if ( ! dat.isValid() ) {
139  cout << myname << "ERROR: Unable to initialize output data container." << endl;
140  return 4;
141  }
142  float fndat = ndatIn;
143  float nfac = dft.normalization().isStandard() ? 1.0/fndat :
144  dft.normalization().isConsistent() ? 1.0/sqrt(fndat) :
145  dft.normalization().isBin() ? 1.0 : 0.0;
146  Plan& plan = backwardPlan(nsams);
147  for ( Index idat=0; idat<ndatOut; ++idat ) floatOutData()[idat] = nfac*dft.floatData()[idat];
148  fftw_execute(plan);
149  dat.copyDataIn(m_inData);
150  return 0;
151 }
Index checkDataSize(const IndexArray &nsams) const
Definition: Fw2dFFT.cxx:39
Plan & backwardPlan(const IndexArray &nsams)
Definition: Fw2dFFT.cxx:77
unsigned int Index
virtual bool isValid() const
Definition: RealDftData.h:42
DftFloat * m_inData
Definition: Fw2dFFT.h:92
static Index dataSize(const std::array< Index, N > &nsams)
Definition: Real2dDftData.h:34
static Index dftFloatDataSize(const std::array< Index, N > &nsams)
const Norm & normalization() const override
Index m_ndatMax
Definition: Fw2dFFT.h:90
QTextStream & endl(QTextStream &s)
fftw_plan Plan
Definition: Fw2dFFT.h:41
DftFloat * floatOutData()
Definition: Fw2dFFT.h:86
int Fw2dFFT::fftForward ( const Data dat,
DFT dft,
Index  logLevel = 0 
)

Definition at line 93 of file Fw2dFFT.cxx.

93  {
94  const string myname = "Fw2dFFT::fftForward: ";
95  if ( ! dat.isValid() ) return 1;
96  auto nsams = dat.nSamples();
97  dft.reset(nsams);
98  if ( ! dft.isValid() ) return 2;
99  if ( dft.normalization().isPower() ) {
100  cout << myname << "ERROR: Power normalization is not (yet) supported." << endl;
101  return 3;
102  }
103  Index ndatIn = dat.size();
104  Index ndatOut = DFT::dftFloatDataSize(nsams);
105  Plan& plan = forwardPlan(nsams);
106  float fndat = ndatIn;
107  float nfac = dft.normalization().isStandard() ? 1.0 :
108  dft.normalization().isConsistent() ? 1.0/sqrt(fndat) :
109  dft.normalization().isBin() ? 1.0/fndat : 0.0;
110  if ( dat.copyDataOut(m_inData) ) {
111  cout << myname << "ERROR: Copy of input data failed." << endl;
112  return 4;
113  }
114  fftw_execute(plan);
115  dft.reset(nsams);
116  for ( Index idat=0; idat<ndatOut; ++idat ) dft.floatData()[idat] = nfac*floatOutData()[idat];
117  return 0;
118 }
auto size() const
unsigned int Index
void reset(Index nsam) override
virtual bool isValid() const
Definition: RealDftData.h:42
DftFloat * m_inData
Definition: Fw2dFFT.h:92
Plan & forwardPlan(const IndexArray &nsams)
Definition: Fw2dFFT.cxx:62
static Index dftFloatDataSize(const std::array< Index, N > &nsams)
const Norm & normalization() const override
QTextStream & endl(QTextStream &s)
fftw_plan Plan
Definition: Fw2dFFT.h:41
DftFloat * floatOutData()
Definition: Fw2dFFT.h:86
FftwComplex* Fw2dFFT::fftwOutData ( )
inline

Definition at line 83 of file Fw2dFFT.h.

83 { return reinterpret_cast<FftwComplex*>(m_outData); }
double FftwComplex[2]
Definition: Fw2dFFT.h:43
Complex * m_outData
Definition: Fw2dFFT.h:93
DftFloat* Fw2dFFT::floatOutData ( )
inline

Definition at line 86 of file Fw2dFFT.h.

86 { return reinterpret_cast<DftFloat*>(m_outData); }
Complex * m_outData
Definition: Fw2dFFT.h:93
double DftFloat
Definition: Fw2dFFT.h:35
Fw2dFFT::Plan & Fw2dFFT::forwardPlan ( const IndexArray nsams)

Definition at line 62 of file Fw2dFFT.cxx.

62  {
63  const string myname = "Fw2dFFT::forwardPlan: ";
64  static Plan badplan;
65  if ( checkDataSize(nsams) ) {
66  cout << myname << "Data cannot be accomodated. Maximum data size is " << m_ndatMax << endl;
67  return badplan;
68  }
69  if ( m_forwardPlans.count(nsams) == 0 ) {
70  m_forwardPlans[nsams] = fftw_plan_dft_r2c_2d(nsams[0], nsams[1], m_inData, fftwOutData(), m_flag);
71  }
72  return m_forwardPlans[nsams];
73 }
Index checkDataSize(const IndexArray &nsams) const
Definition: Fw2dFFT.cxx:39
Index m_flag
Definition: Fw2dFFT.h:91
PlanMap m_forwardPlans
Definition: Fw2dFFT.h:94
FftwComplex * fftwOutData()
Definition: Fw2dFFT.h:83
DftFloat * m_inData
Definition: Fw2dFFT.h:92
Index m_ndatMax
Definition: Fw2dFFT.h:90
QTextStream & endl(QTextStream &s)
fftw_plan Plan
Definition: Fw2dFFT.h:41
bool Fw2dFFT::haveBackwardPlan ( const IndexArray nsams) const

Definition at line 56 of file Fw2dFFT.cxx.

56  {
57  return m_backwardPlans.count(nsams);
58 }
PlanMap m_backwardPlans
Definition: Fw2dFFT.h:95
bool Fw2dFFT::haveForwardPlan ( const IndexArray nsams) const

Definition at line 50 of file Fw2dFFT.cxx.

50  {
51  return m_forwardPlans.count(nsams);
52 }
PlanMap m_forwardPlans
Definition: Fw2dFFT.h:94

Member Data Documentation

PlanMap Fw2dFFT::m_backwardPlans
private

Definition at line 95 of file Fw2dFFT.h.

Index Fw2dFFT::m_flag
private

Definition at line 91 of file Fw2dFFT.h.

PlanMap Fw2dFFT::m_forwardPlans
private

Definition at line 94 of file Fw2dFFT.h.

DftFloat* Fw2dFFT::m_inData
private

Definition at line 92 of file Fw2dFFT.h.

Index Fw2dFFT::m_ndatMax
private

Definition at line 90 of file Fw2dFFT.h.

Complex* Fw2dFFT::m_outData
private

Definition at line 93 of file Fw2dFFT.h.


The documentation for this class was generated from the following files: