Page MenuHomec4science

fftransform_fftw.cpp
No OneTemporary

File Metadata

Created
Tue, Jun 25, 16:44

fftransform_fftw.cpp

/**
* @file
* @section LICENSE
*
* Copyright (©) 2016-19 EPFL (École Polytechnique Fédérale de Lausanne),
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Affero General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU Affero General Public License for more details.
*
* You should have received a copy of the GNU Affero General Public License
* along with this program. If not, see <https://www.gnu.org/licenses/>.
*
*/
/* -------------------------------------------------------------------------- */
#include "fftransform_fftw.hh"
#define TAMAAS_FFTW_FLAG FFTW_ESTIMATE
namespace tamaas {
template <typename T, UInt dim>
FFTransformFFTW<T, dim>::FFTransformFFTW(Grid<T, dim>& real,
GridHermitian<T, dim>& spectral)
: FFTransform<T, dim>(real, spectral), forward_plan(nullptr),
backward_plan(nullptr) {
TAMAAS_ASSERT(real.getNbComponents() == spectral.getNbComponents(),
"components for FFTW do not match");
// Data pointers
Real* in = this->real.getInternalData();
complex_t* out =
reinterpret_cast<complex_t*>(this->spectral.getInternalData());
int components = real.getNbComponents();
// fftw parameters
int rank = dim; // dimension of fft
int n[dim] = {0};
std::copy(real.sizes().begin(), real.sizes().end(),
n); // size of individual fft
int howmany = components; // how many fft to compute
int idist = 1, odist = 1; // components are next to each other in memory
int istride = real.getStrides().back() * components,
ostride = spectral.getStrides().back() * components;
int *inembed = nullptr, *onembed = nullptr; // row major
#if TAMAAS_FFTW_FLAG != FFTW_ESTIMATE
// backup the input
Grid<Real, dim> backup(real);
#endif
forward_plan =
fftw::plan_many_forward(rank, n, howmany, in, inembed, istride, idist,
out, onembed, ostride, odist, TAMAAS_FFTW_FLAG);
backward_plan =
fftw::plan_many_backward(rank, n, howmany, out, onembed, ostride, odist,
in, inembed, istride, idist, TAMAAS_FFTW_FLAG);
#if TAMAAS_FFTW_FLAG != FFTW_ESTIMATE
// restore the backup
real = backup;
#endif
}
/* -------------------------------------------------------------------------- */
template <typename T, UInt dim>
FFTransformFFTW<T, dim>::~FFTransformFFTW() {
fftw::destroy(forward_plan);
fftw::destroy(backward_plan);
}
/* -------------------------------------------------------------------------- */
template <typename T, UInt dim>
void FFTransformFFTW<T, dim>::forwardTransform() {
fftw::execute(forward_plan);
}
/* -------------------------------------------------------------------------- */
template <typename T, UInt dim>
void FFTransformFFTW<T, dim>::backwardTransform() {
fftw::execute(backward_plan);
this->normalize();
}
/* -------------------------------------------------------------------------- */
template class FFTransformFFTW<Real, 1>;
template class FFTransformFFTW<Real, 2>;
template class FFTransformFFTW<Real, 3>;
} // namespace tamaas

Event Timeline