LSSTApplications  10.0+286,10.0+36,10.0+46,10.0-2-g4f67435,10.1+152,10.1+37,11.0,11.0+1,11.0-1-g47edd16,11.0-1-g60db491,11.0-1-g7418c06,11.0-2-g04d2804,11.0-2-g68503cd,11.0-2-g818369d,11.0-2-gb8b8ce7
LSSTDataManagementBasePackage
Public Types | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes | List of all members
ndarray::FourierTransform< T, N > Class Template Reference

A wrapper for FFTW plans for fast Fourier transforms. More...

#include <FourierTransform.h>

Inheritance diagram for ndarray::FourierTransform< T, N >:

Public Types

typedef boost::shared_ptr
< FourierTransform
Ptr
 
typedef detail::FourierTraits
< T >::ElementX 
ElementX
 Real-space array data type;. More...
 
typedef detail::FourierTraits
< T >::ElementK 
ElementK
 Fourier-space array data type;. More...
 
typedef Vector< int, N > Index
 Shape type for arrays. More...
 
typedef Array< ElementX, N, N > ArrayX
 Real-space array type. More...
 
typedef Array< ElementK, N, N > ArrayK
 Fourier-space array type. More...
 
typedef Vector< int, N+1 > MultiplexIndex
 Shape type for multiplexed arrays. More...
 
typedef Array< ElementX, N+1, N+1 > MultiplexArrayX
 Real-space multiplexed array type. More...
 
typedef Array< ElementK, N+1, N+1 > MultiplexArrayK
 Fourier-space multiplexed array type. More...
 

Public Member Functions

void execute ()
 Execute the FFTW plan. More...
 
 ~FourierTransform ()
 
template<int M>
Array< typename
FourierTransform< T, N >
::ElementX, M, M > 
initializeX (Vector< int, M > const &shape)
 
template<int M>
Array< typename
FourierTransform< T, N >
::ElementK, M, M > 
initializeK (Vector< int, M > const &shape)
 

Static Public Member Functions

static Ptr planForward (Index const &shape, ArrayX &x, ArrayK &k)
 Create a plan for forward-transforming a single N-dimensional array. More...
 
static Ptr planInverse (Index const &shape, ArrayK &k, ArrayX &x)
 Create a plan for inverse-transforming a single N-dimensional array. More...
 
static Ptr planMultiplexForward (MultiplexIndex const &shape, MultiplexArrayX &x, MultiplexArrayK &k)
 Create a plan for forward-transforming a sequence of nested N-dimensional arrays. More...
 
static Ptr planMultiplexInverse (MultiplexIndex const &shape, MultiplexArrayK &k, MultiplexArrayX &x)
 Create a plan for inverse-transforming a sequence of nested N-dimensional arrays. More...
 
template<int M>
static Array< ElementX, M, M > initializeX (Vector< int, M > const &shape)
 Create a new real-space array with the given real-space shape. More...
 
template<int M>
static Array< ElementK, M, M > initializeK (Vector< int, M > const &shape)
 Create a new Fourier-space array with the given real-space shape. More...
 
template<int M>
static void initialize (Vector< int, M > const &shape, Array< ElementX, M, M > &x, Array< ElementK, M, M > &k)
 Initialize, as necessary, a pair of arrays with the given real-space shape. More...
 

Private Types

typedef boost::shared_ptr
< ElementX
OwnerX
 
typedef boost::shared_ptr
< ElementK
OwnerK
 

Private Member Functions

 BOOST_STATIC_ASSERT ((!boost::is_const< T >::value))
 
 FourierTransform (void *plan, Manager::Ptr const &x, Manager::Ptr const &k)
 

Private Attributes

void * _plan
 
Manager::Ptr _x
 
Manager::Ptr _k
 

Detailed Description

template<typename T, int N>
class ndarray::FourierTransform< T, N >

A wrapper for FFTW plans for fast Fourier transforms.

An instance of FourierTransform holds an FFTW "plan", providing repeated forward or inverse FFTs of predetermined arrays.

Multiplex plans can also be generated to perform an N-dimensional FFT on the nested arrays of an (N+1)-dimensional array.

Static member functions of FourierTransform are used to create instances, and optionally initialize the involved arrays.

Definition at line 54 of file FourierTransform.h.

Member Typedef Documentation

template<typename T, int N>
typedef Array<ElementK,N,N> ndarray::FourierTransform< T, N >::ArrayK

Fourier-space array type.

Definition at line 65 of file FourierTransform.h.

template<typename T, int N>
typedef Array<ElementX,N,N> ndarray::FourierTransform< T, N >::ArrayX

Real-space array type.

Definition at line 64 of file FourierTransform.h.

template<typename T, int N>
typedef detail::FourierTraits<T>::ElementK ndarray::FourierTransform< T, N >::ElementK

Fourier-space array data type;.

Definition at line 61 of file FourierTransform.h.

template<typename T, int N>
typedef detail::FourierTraits<T>::ElementX ndarray::FourierTransform< T, N >::ElementX

Real-space array data type;.

Definition at line 60 of file FourierTransform.h.

template<typename T, int N>
typedef Vector<int,N> ndarray::FourierTransform< T, N >::Index

Shape type for arrays.

Definition at line 63 of file FourierTransform.h.

template<typename T, int N>
typedef Array<ElementK,N+1,N+1> ndarray::FourierTransform< T, N >::MultiplexArrayK

Fourier-space multiplexed array type.

Definition at line 68 of file FourierTransform.h.

template<typename T, int N>
typedef Array<ElementX,N+1,N+1> ndarray::FourierTransform< T, N >::MultiplexArrayX

Real-space multiplexed array type.

Definition at line 67 of file FourierTransform.h.

template<typename T, int N>
typedef Vector<int,N+1> ndarray::FourierTransform< T, N >::MultiplexIndex

Shape type for multiplexed arrays.

Definition at line 66 of file FourierTransform.h.

template<typename T, int N>
typedef boost::shared_ptr<ElementK> ndarray::FourierTransform< T, N >::OwnerK
private

Definition at line 141 of file FourierTransform.h.

template<typename T, int N>
typedef boost::shared_ptr<ElementX> ndarray::FourierTransform< T, N >::OwnerX
private

Definition at line 140 of file FourierTransform.h.

template<typename T, int N>
typedef boost::shared_ptr<FourierTransform> ndarray::FourierTransform< T, N >::Ptr

Definition at line 58 of file FourierTransform.h.

Constructor & Destructor Documentation

template<typename T , int N>
ndarray::FourierTransform< T, N >::~FourierTransform ( )

Definition at line 156 of file FourierTransform.cc.

156  {
157  detail::FFTWTraits<T>::destroy(
158  reinterpret_cast<typename detail::FFTWTraits<T>::Plan>(_plan)
159  );
160 }
template<typename T, int N>
ndarray::FourierTransform< T, N >::FourierTransform ( void *  plan,
Manager::Ptr const &  x,
Manager::Ptr const &  k 
)
inlineprivate

Definition at line 143 of file FourierTransform.h.

144  : _plan(plan), _x(x), _k(k) {}
double x

Member Function Documentation

template<typename T, int N>
ndarray::FourierTransform< T, N >::BOOST_STATIC_ASSERT ( (!boost::is_const< T >::value)  )
private
template<typename T , int N>
void ndarray::FourierTransform< T, N >::execute ( )

Execute the FFTW plan.

Definition at line 149 of file FourierTransform.cc.

149  {
150  detail::FFTWTraits<T>::execute(
151  reinterpret_cast<typename detail::FFTWTraits<T>::Plan>(_plan)
152  );
153 }
template<typename T , int N>
template<int M>
void ndarray::FourierTransform< T, N >::initialize ( Vector< int, M > const &  shape,
Array< ElementX, M, M > &  x,
Array< ElementK, M, M > &  k 
)
static

Initialize, as necessary, a pair of arrays with the given real-space shape.

If either array is not empty, it must be consistent with the given shape.

Definition at line 49 of file FourierTransform.cc.

53  {
54  if (x.empty()) x = initializeX(shape);
55  if (k.empty()) k = initializeK(shape);
56  NDARRAY_ASSERT(x.getShape() == shape);
57  NDARRAY_ASSERT(std::equal(shape.begin(), shape.end()-1, k.getShape().begin()));
58 }
static Array< ElementX, M, M > initializeX(Vector< int, M > const &shape)
Create a new real-space array with the given real-space shape.
#define NDARRAY_ASSERT(ARG)
Definition: ndarray_fwd.h:51
static Array< ElementK, M, M > initializeK(Vector< int, M > const &shape)
Create a new Fourier-space array with the given real-space shape.
double x
< unspecified-expression-type > equal(ExpressionBase< Operand > const &operand, Scalar const &scalar)
Definition: operators.h:796
template<typename T, int N>
template<int M>
Array<typename FourierTransform<T,N>::ElementK,M,M> ndarray::FourierTransform< T, N >::initializeK ( Vector< int, M > const &  shape)

Definition at line 39 of file FourierTransform.cc.

39  {
40  Vector<int,M> kShape(shape);
41  kShape[M-1] = detail::FourierTraits<T>::computeLastDimensionSize(shape[M-1]);
42  OwnerK kOwner = detail::FFTWTraits<T>::allocateK(kShape.product());
43  return Array<ElementK,M,M>(external(kOwner.get(), kShape, ROW_MAJOR, kOwner));
44 }
detail::ExternalInitializer< T, N, Owner > external(T *data, Vector< int, N > const &shape, Vector< int, N > const &strides, Owner const &owner)
Create an expression that initializes an Array with externally allocated memory.
boost::shared_ptr< ElementK > OwnerK
template<typename T, int N>
template<int M>
static Array<ElementK,M,M> ndarray::FourierTransform< T, N >::initializeK ( Vector< int, M > const &  shape)
static

Create a new Fourier-space array with the given real-space shape.

template<typename T, int N>
template<int M>
Array<typename FourierTransform<T,N>::ElementX,M,M> ndarray::FourierTransform< T, N >::initializeX ( Vector< int, M > const &  shape)

Definition at line 31 of file FourierTransform.cc.

31  {
32  OwnerX xOwner = detail::FFTWTraits<T>::allocateX(shape.product());
33  return Array<ElementX,M,M>(external(xOwner.get(), shape, ROW_MAJOR, xOwner));
34 }
boost::shared_ptr< ElementX > OwnerX
detail::ExternalInitializer< T, N, Owner > external(T *data, Vector< int, N > const &shape, Vector< int, N > const &strides, Owner const &owner)
Create an expression that initializes an Array with externally allocated memory.
template<typename T, int N>
template<int M>
static Array<ElementX,M,M> ndarray::FourierTransform< T, N >::initializeX ( Vector< int, M > const &  shape)
static

Create a new real-space array with the given real-space shape.

template<typename T, int N>
FourierTransform< T, N >::Ptr ndarray::FourierTransform< T, N >::planForward ( Index const &  shape,
ArrayX x,
ArrayK k 
)
static

Create a plan for forward-transforming a single N-dimensional array.

Arrays will be initialized with new memory if empty. If they are not empty, existing data may be overwritten when the plan is created.

Parameters
shapeShape of the real-space array.
xInput real-space array.
kOutput Fourier-space array.

Definition at line 62 of file FourierTransform.cc.

66  {
67  initialize(shape,x,k);
68  return Ptr(
69  new FourierTransform(
70  detail::FFTWTraits<T>::forward(
71  N, shape.begin(), 1,
72  x.getData(), NULL, 1, 0,
73  k.getData(), NULL, 1, 0,
74  FFTW_MEASURE | FFTW_DESTROY_INPUT
75  ),
76  x.getManager(),
77  k.getManager()
78  )
79  );
80 }
static void initialize(Vector< int, M > const &shape, Array< ElementX, M, M > &x, Array< ElementK, M, M > &k)
Initialize, as necessary, a pair of arrays with the given real-space shape.
FourierTransform(void *plan, Manager::Ptr const &x, Manager::Ptr const &k)
double x
boost::shared_ptr< FourierTransform > Ptr
template<typename T, int N>
FourierTransform< T, N >::Ptr ndarray::FourierTransform< T, N >::planInverse ( Index const &  shape,
ArrayK k,
ArrayX x 
)
static

Create a plan for inverse-transforming a single N-dimensional array.

Arrays will be initialized with new memory if empty. If they are not empty, existing data may be overwritten when the plan is created.

Parameters
shapeShape of the real-space array.
kInput Fourier-space array.
xOutput real-space array.

Definition at line 84 of file FourierTransform.cc.

88  {
89  initialize(shape,x,k);
90  return Ptr(
91  new FourierTransform(
92  detail::FFTWTraits<T>::inverse(
93  N, shape.begin(), 1,
94  k.getData(), NULL, 1, 0,
95  x.getData(), NULL, 1, 0,
96  FFTW_MEASURE | FFTW_DESTROY_INPUT
97  ),
98  x.getManager(),
99  k.getManager()
100  )
101  );
102 }
static void initialize(Vector< int, M > const &shape, Array< ElementX, M, M > &x, Array< ElementK, M, M > &k)
Initialize, as necessary, a pair of arrays with the given real-space shape.
FourierTransform(void *plan, Manager::Ptr const &x, Manager::Ptr const &k)
double x
boost::shared_ptr< FourierTransform > Ptr
template<typename T, int N>
FourierTransform< T, N >::Ptr ndarray::FourierTransform< T, N >::planMultiplexForward ( MultiplexIndex const &  shape,
MultiplexArrayX x,
MultiplexArrayK k 
)
static

Create a plan for forward-transforming a sequence of nested N-dimensional arrays.

Arrays will be initialized with new memory if empty. If they are not empty, existing data may be overwritten when the plan is created.

Parameters
shapeShape of the real-space array. First dimension is multiplexed.
xInput real-space array.
kOutput Fourier-space array.

Definition at line 106 of file FourierTransform.cc.

110  {
111  initialize(shape,x,k);
112  return Ptr(
113  new FourierTransform(
114  detail::FFTWTraits<T>::forward(
115  N, shape.begin()+1, shape[0],
116  x.getData(), NULL, 1, x.template getStride<0>(),
117  k.getData(), NULL, 1, k.template getStride<0>(),
118  FFTW_MEASURE | FFTW_DESTROY_INPUT
119  ),
120  x.getManager(),
121  k.getManager()
122  )
123  );
124 }
static void initialize(Vector< int, M > const &shape, Array< ElementX, M, M > &x, Array< ElementK, M, M > &k)
Initialize, as necessary, a pair of arrays with the given real-space shape.
FourierTransform(void *plan, Manager::Ptr const &x, Manager::Ptr const &k)
double x
boost::shared_ptr< FourierTransform > Ptr
template<typename T, int N>
FourierTransform< T, N >::Ptr ndarray::FourierTransform< T, N >::planMultiplexInverse ( MultiplexIndex const &  shape,
MultiplexArrayK k,
MultiplexArrayX x 
)
static

Create a plan for inverse-transforming a sequence of nested N-dimensional arrays.

Arrays will be initialized with new memory if empty. If they are not empty, existing data may be overwritten when the plan is created.

Parameters
shapeShape of the real-space array. First dimension is multiplexed.
kInput Fourier-space array.
xOutput real-space array.

Definition at line 128 of file FourierTransform.cc.

132  {
133  initialize(shape,x,k);
134  return Ptr(
135  new FourierTransform(
136  detail::FFTWTraits<T>::inverse(
137  N, shape.begin()+1, shape[0],
138  k.getData(), NULL, 1, k.template getStride<0>(),
139  x.getData(), NULL, 1, x.template getStride<0>(),
140  FFTW_MEASURE | FFTW_DESTROY_INPUT
141  ),
142  x.getManager(),
143  k.getManager()
144  )
145  );
146 }
static void initialize(Vector< int, M > const &shape, Array< ElementX, M, M > &x, Array< ElementK, M, M > &k)
Initialize, as necessary, a pair of arrays with the given real-space shape.
FourierTransform(void *plan, Manager::Ptr const &x, Manager::Ptr const &k)
double x
boost::shared_ptr< FourierTransform > Ptr

Member Data Documentation

template<typename T, int N>
Manager::Ptr ndarray::FourierTransform< T, N >::_k
private

Definition at line 148 of file FourierTransform.h.

template<typename T, int N>
void* ndarray::FourierTransform< T, N >::_plan
private

Definition at line 146 of file FourierTransform.h.

template<typename T, int N>
Manager::Ptr ndarray::FourierTransform< T, N >::_x
private

Definition at line 147 of file FourierTransform.h.


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