CppNoddy  0.92
Loading...
Searching...
No Matches
Functions
CppNoddy::FT Namespace Reference

Some algorithms associated with Fourier transforms of CppNoddy container data. More...

Functions

OneD_Node_Mesh< D_complexshift (const OneD_Node_Mesh< D_complex > &ft)
 Shift the frequency spectrum obtained from dft to give positive and negative freq. More...
 
OneD_Node_Mesh< D_complexishift (const OneD_Node_Mesh< D_complex > &ft)
 Invert the shif operation to recover a spectrum in the form that is expected by the idft routine. More...
 
OneD_Node_Mesh< D_complexdft_with_shift (const OneD_Node_Mesh< D_complex > &f)
 A wrapper that calls the 'dft' method above followed by the 'shift' method. More...
 
TwoD_Node_Mesh< D_complexdft_with_shift (const TwoD_Node_Mesh< D_complex > &f)
 
OneD_Node_Mesh< D_complexidft_with_ishift (const OneD_Node_Mesh< D_complex > &ft, double origin=0)
 A wrapper that calls the 'ishift' method above followed by the 'idft' method. More...
 
TwoD_Node_Mesh< D_complexidft_with_ishift (const TwoD_Node_Mesh< D_complex > &ft, double origin=0)
 
OneD_Node_Mesh< D_complexdft (const OneD_Node_Mesh< D_complex > &f)
 (Slow) DFT of the real data (x_i,f_i), i = 0, ... N-1; N must be EVEN. More...
 
OneD_Node_Mesh< D_complexidft (const OneD_Node_Mesh< D_complex > &ft, double origin=0)
 (Slow) Inverse DFT of the complex data (omega_i,F_i), i = 0,...,N-1; N must be EVEN. More...
 

Detailed Description

Some algorithms associated with Fourier transforms of CppNoddy container data.

Function Documentation

◆ dft()

OneD_Node_Mesh< D_complex > CppNoddy::FT::dft ( const OneD_Node_Mesh< D_complex > &  f)

(Slow) DFT of the real data (x_i,f_i), i = 0, ... N-1; N must be EVEN.

Parameters
fThe data to be (discrete) Fourier transformed.
Returns
The Fourier transform in "frequency" space (omega_i,F_i). The N frequencies are omega_i = 0, domega,..., (N-1)*domega. The SECOND half of the spectrum provides the -ve frequencies (as in MATLAB).

Definition at line 88 of file FT.cpp.

88 {
89 const D_complex I(0.0,1.0); // imaginary unit
90 // number of samples in the signal
91 std::size_t N = f.get_nnodes();
92 // time step in the signal (has to be constant)
93 double dt = f.coord(1)-f.coord(0);
94 // sample frequency (rad) is related to inverse of "time" step
95 double fs = 2*M_PI/dt;
96 // frequency step between frequency nodes
97 // - largest frequency is (N-1)*df=(N-1)*fs/N
98 double df = fs/N;
99 // container to return the Fourier spectrum
100 OneD_Node_Mesh<D_complex> ft( f.nodes(), f.get_nvars() );
101 //
102 const double Wn = 2*M_PI/N;
103 // slow FT
104 for ( std::size_t k = 0; k < N; ++k ) {
105 // set the frequency
106 ft.coord(k) = k*df; // 0,df,...,(N-1)*df
107 //
108 //
109 // t = n*dt
110 // omega = k*df
111 // df = (2*pi/dt)/N
112 //
113 // variables are all initialised to zero
114 DenseVector<D_complex> vars( f.get_nvars(), 0.0 );
115 for ( std::size_t n = 0; n < N; ++n ) {
116 vars += f.get_nodes_vars(n)*exp(-I*Wn*double(k*n));
117 }
118 // multiply by "dt" to make this correlate with analytical
119 // results from the FT integral
120 ft.set_nodes_vars(k,vars*dt);
121 }
122 return ft;
123 }
@ f
Definition: BVPBerman.cpp:15
An DenseVector class – a dense vector object.
Definition: DenseVector.h:34
A one dimensional mesh utility object.
std::complex< double > D_complex
A complex double precision number using std::complex.
Definition: Types.h:98

References CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::coord(), f, and CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::set_nodes_vars().

Referenced by dft_with_shift(), and main().

◆ dft_with_shift() [1/2]

OneD_Node_Mesh< D_complex > CppNoddy::FT::dft_with_shift ( const OneD_Node_Mesh< D_complex > &  f)

A wrapper that calls the 'dft' method above followed by the 'shift' method.

Returns
The Fourier transform in "frequency" space (omega_i,F_i). The N frequencies are omega_i = -(N/2)*domega, ..., -domega , 0 , +domega, ... ,+(N/2-1)*domega.

Definition at line 48 of file FT.cpp.

48 {
49 return shift(dft(f));
50 }
OneD_Node_Mesh< D_complex > dft(const OneD_Node_Mesh< D_complex > &f)
(Slow) DFT of the real data (x_i,f_i), i = 0, ... N-1; N must be EVEN.
Definition: FT.cpp:88
OneD_Node_Mesh< D_complex > shift(const OneD_Node_Mesh< D_complex > &ft)
Shift the frequency spectrum obtained from dft to give positive and negative freq.
Definition: FT.cpp:15

References dft(), f, and shift().

Referenced by dft_with_shift(), and main().

◆ dft_with_shift() [2/2]

TwoD_Node_Mesh< D_complex > CppNoddy::FT::dft_with_shift ( const TwoD_Node_Mesh< D_complex > &  f)

Definition at line 52 of file FT.cpp.

52 {
53 // assume that we idft each ft(.,node)
54 unsigned ny = f.get_nnodes().second;
55 unsigned nx = f.get_nnodes().first;
56 TwoD_Node_Mesh<D_complex> FTsoln( f.xnodes(), f.ynodes(), f.get_nvars());
57 for ( unsigned jy = 0; jy < ny; ++jy ) {
58 OneD_Node_Mesh<D_complex> mesh = f.get_xsection_at_ynode(jy);
59 OneD_Node_Mesh<D_complex> ftmesh = FT::dft_with_shift(mesh);
60 for ( unsigned jx = 0; jx < nx; ++jx ) {
61 FTsoln.xcoord(jx)=ftmesh.coord(jx);
62 FTsoln.set_nodes_vars( jx, jy, ftmesh.get_nodes_vars(jx) );
63 }
64 }
65 return FTsoln;
66 }
DenseVector< _Type > get_nodes_vars(const std::size_t &node) const
Get the variables stored at A SPECIFIED node.
const _Xtype & coord(const std::size_t &node) const
Access a nodal position.
A two dimensional mesh utility object.

References CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::coord(), dft_with_shift(), f, CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::get_nodes_vars(), CppNoddy::TwoD_Node_Mesh< _Type >::set_nodes_vars(), and CppNoddy::TwoD_Node_Mesh< _Type >::xcoord().

◆ idft()

OneD_Node_Mesh< D_complex > CppNoddy::FT::idft ( const OneD_Node_Mesh< D_complex > &  ft,
double  origin = 0 
)

(Slow) Inverse DFT of the complex data (omega_i,F_i), i = 0,...,N-1; N must be EVEN.

Parameters
ftThe frequency spectrum in the (unshifted) format.
originThe position of the first data point returned.
Returns
The inversion data (origin+x_i, f_i).

Definition at line 153 of file FT.cpp.

153 {
154 const D_complex I(0.0,1.0); // imaginary unit
155 // number of "frequencies" in the data
156 std::size_t N = ft.get_nnodes();
157 // (uniform) frequency step in the data
158 double df = ft.coord(1)-ft.coord(0);
159 // sample frequency
160 double fs = N*df;
161 // time step in output signal
162 double dt = 2*M_PI/fs;
163 // NOTE: dt*df = 2*pi/N
164 //
165 // containter to return the "time series" in
167 //
168 const double Wn = 2*M_PI/N;
169 // number of variables at each time/frequency
170 const std::size_t nVars( ft.get_nvars() );
171 // slow inverse FT
172 for ( std::size_t n = 0; n < N; ++n ) {
173 // set the coordinate/time using the offset passed to this method
174 f.coord(n) = origin + n*dt;
175 // invert for all variables in the mesh at the same time
176 DenseVector<D_complex> vars( nVars, 0.0 );
177 for ( std::size_t k = 0; k < N; ++k ) {
178 vars += ft.get_nodes_vars(k)*exp(+I*Wn*double(k*n));
179 }
180 // the dt factor is for consistency with continuous FT integral
181 f.set_nodes_vars(n,vars/(N*dt));
182 }
183 return f;
184 }
std::size_t get_nvars() const
std::size_t get_nnodes() const

References CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::coord(), f, CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::get_nnodes(), CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::get_nodes_vars(), and CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::get_nvars().

Referenced by idft_with_ishift(), and main().

◆ idft_with_ishift() [1/2]

OneD_Node_Mesh< D_complex > CppNoddy::FT::idft_with_ishift ( const OneD_Node_Mesh< D_complex > &  ft,
double  origin = 0 
)

A wrapper that calls the 'ishift' method above followed by the 'idft' method.

Parameters
ftThe frequency spectrum in the (unshifted) format.
originThe position of the first data point returned.
Returns
The inversion data (origin+x_i, f_i).

Definition at line 68 of file FT.cpp.

68 {
69 return idft(ishift(ft),origin);
70 }
OneD_Node_Mesh< D_complex > ishift(const OneD_Node_Mesh< D_complex > &ft)
Invert the shif operation to recover a spectrum in the form that is expected by the idft routine.
Definition: FT.cpp:31
OneD_Node_Mesh< D_complex > idft(const OneD_Node_Mesh< D_complex > &ft, double origin=0)
(Slow) Inverse DFT of the complex data (omega_i,F_i), i = 0,...,N-1; N must be EVEN.
Definition: FT.cpp:153

References idft(), and ishift().

Referenced by idft_with_ishift(), and main().

◆ idft_with_ishift() [2/2]

TwoD_Node_Mesh< D_complex > CppNoddy::FT::idft_with_ishift ( const TwoD_Node_Mesh< D_complex > &  ft,
double  origin = 0 
)

Definition at line 72 of file FT.cpp.

72 {
73 // assume that we idft each ft(.,node)
74 unsigned ny = ft.get_nnodes().second;
75 unsigned nx = ft.get_nnodes().first;
76 TwoD_Node_Mesh<D_complex> soln( ft.xnodes(), ft.ynodes(), ft.get_nvars());
77 for ( unsigned jy = 0; jy < ny; ++jy ) {
79 OneD_Node_Mesh<D_complex> mesh = FT::idft_with_ishift(ftmesh,origin);
80 for ( unsigned jx = 0; jx < nx; ++jx ) {
81 soln.xcoord(jx)=mesh.coord(jx);
82 soln.set_nodes_vars( jx, jy, mesh.get_nodes_vars(jx) );
83 }
84 }
85 return soln;
86 }
std::pair< std::size_t, std::size_t > get_nnodes() const
Get the number of nodes in the two directions of the 2D mesh.
const DenseVector< double > & xnodes() const
Access the vector of x-nodal positions.
std::size_t get_nvars() const
Get the number of variables that are stored at each node.
const DenseVector< double > & ynodes() const
Access the vector of y-nodal positions.
OneD_Node_Mesh< _Type > get_xsection_at_ynode(const std::size_t nodey) const
Get a cross section of the 2D mesh at a specified (constant) y node.

References CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::coord(), CppNoddy::TwoD_Node_Mesh< _Type >::get_nnodes(), CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::get_nodes_vars(), CppNoddy::TwoD_Node_Mesh< _Type >::get_nvars(), CppNoddy::TwoD_Node_Mesh< _Type >::get_xsection_at_ynode(), idft_with_ishift(), CppNoddy::TwoD_Node_Mesh< _Type >::set_nodes_vars(), CppNoddy::TwoD_Node_Mesh< _Type >::xcoord(), CppNoddy::TwoD_Node_Mesh< _Type >::xnodes(), and CppNoddy::TwoD_Node_Mesh< _Type >::ynodes().

◆ ishift()

OneD_Node_Mesh< D_complex > CppNoddy::FT::ishift ( const OneD_Node_Mesh< D_complex > &  ft)

Invert the shif operation to recover a spectrum in the form that is expected by the idft routine.

Parameters
ftA previously shifted frequency spectrum.
Returns
The frequency spectrum in the (unshifted) format. 0,domega,...,(N/2)*domega, (N/2+1)*domega,...,(N-1)*domega.

Definition at line 31 of file FT.cpp.

31 {
32 OneD_Node_Mesh<D_complex> ft_ishifted( ft );
33 double df = ft.coord(1)-ft.coord(0);
34 std::size_t N = ft.get_nnodes();
35 for ( std::size_t i = 0; i < N/2; ++i ) {
36 for ( std::size_t var = 0; var < ft.get_nvars(); ++var ) {
37 // first half of shifted output
38 ft_ishifted( i, var ) = ft( N/2 + i, var );
39 // second half of shifted output
40 ft_ishifted( N/2 + i, var ) = ft( i, var );
41 }
42 ft_ishifted.coord( i ) = ft.coord( N/2 + i );
43 ft_ishifted.coord( N/2 + i ) = ft.coord(N-1) + i*df;
44 }
45 return ft_ishifted;
46 }

References CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::coord(), CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::get_nnodes(), and CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::get_nvars().

Referenced by idft_with_ishift().

◆ shift()

OneD_Node_Mesh< D_complex > CppNoddy::FT::shift ( const OneD_Node_Mesh< D_complex > &  ft)

Shift the frequency spectrum obtained from dft to give positive and negative freq.

Note: input should have an even number of nodes.

Parameters
ftThe frequency spectrum in the (unshifted) format.
Returns
The frequency spectrum in the (shifted) format. -(N/2)domega,...,-domega,0,domega,...,(N/2-1)*domega

Definition at line 15 of file FT.cpp.

15 {
16 OneD_Node_Mesh<D_complex> ft_shifted( ft );
17 std::size_t N = ft.get_nnodes();
18 for ( std::size_t i = 0; i < N/2; ++i ) {
19 for ( std::size_t var = 0; var < ft.get_nvars(); ++var ) {
20 // first half of shifted output
21 ft_shifted( i, var ) = ft( N/2+i, var );
22 // second half of shifted output
23 ft_shifted( N/2 + i, var ) = ft( i, var );
24 }
25 ft_shifted.coord( i ) = -ft.coord( N/2 - i );
26 ft_shifted.coord( N/2 + i ) = ft.coord( i );
27 }
28 return ft_shifted;
29 }

References CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::coord(), CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::get_nnodes(), and CppNoddy::OneD_Node_Mesh< _Type, _Xtype >::get_nvars().

Referenced by dft_with_shift(), and main().

© 2012

R.E. Hewitt