Skip to content
Snippets Groups Projects
Commit 4eaa71f6 authored by ansari's avatar ansari
Browse files

Ajout documentation - Reza 15/2/2001

parent 00b9a520
No related branches found
No related tags found
No related merge requests found
......@@ -2,9 +2,34 @@
#include "FFTW/fftw.h"
#include "FFTW/rfftw.h"
/*!
\class SOPHYA::FFTWServer
\ingroup IFFTW
An implementation of FFTServerInterface based on FFTW, double
precision arrays, using FFTW package, availabale from http://www.fftw.org.
Refer to FFTServerInterface for details about FFTServer operations.
\code
#include "fftwserver.h"
// ...
TMatrix<r_8> in(24,32);
TMatrix< complex<r_8> > out;
in = RandomSequence();
FFTWServer ffts;
ffts.setNormalize(true); // To have normalized transforms
cout << " FFTServer info string= " << ffts.getInfo() << endl;
cout << "in= " << in << endl;
cout << " Calling ffts.FFTForward(in, out) : " << endl;
ffts.FFTForward(in, out);
cout << "out= " << out << endl;
\endcode
*/
#define MAXND_FFTW 5
namespace SOPHYA {
class FFTWServerPlan {
public:
FFTWServerPlan(int n, fftw_direction dir, bool fgreal=false);
......@@ -27,6 +52,8 @@ public:
};
} // Fin du namespace
FFTWServerPlan::FFTWServerPlan(int n, fftw_direction dir, bool fgreal)
{
if (n < 1)
......@@ -39,8 +66,17 @@ FFTWServerPlan::FFTWServerPlan(int n, fftw_direction dir, bool fgreal)
for(int k=0; k<MAXND_FFTW; k++) _nxyz[k] = -10;
_n = n;
_dir = dir;
if (fgreal) rp = rfftw_create_plan(n, dir, FFTW_ESTIMATE);
else p = fftw_create_plan(n, dir, FFTW_ESTIMATE);
if (fgreal) {
rp = rfftw_create_plan(n, dir, FFTW_ESTIMATE);
if (rp == NULL)
throw AllocationError("FFTWServerPlan: failed to create plan (rp) !");
}
else {
p = fftw_create_plan(n, dir, FFTW_ESTIMATE);
if (p == NULL)
throw AllocationError("FFTWServerPlan: failed to create plan (p) !");
}
}
FFTWServerPlan::FFTWServerPlan(int nd, int * nxyz, fftw_direction dir, bool fgreal)
......@@ -61,8 +97,17 @@ FFTWServerPlan::FFTWServerPlan(int nd, int * nxyz, fftw_direction dir, bool fgre
}
for(k=nd; k<MAXND_FFTW; k++) _nxyz[k] = -10;
_dir = dir;
if (fgreal) rpnd = rfftwnd_create_plan(_nd, _nxyz, dir, FFTW_ESTIMATE);
else pnd = fftwnd_create_plan(_nd, _nxyz, dir, FFTW_ESTIMATE);
if (fgreal) {
rpnd = rfftwnd_create_plan(_nd, _nxyz, dir, FFTW_ESTIMATE);
if (rpnd == NULL)
throw AllocationError("FFTWServerPlan: failed to create plan (rpnd) !");
}
else {
pnd = fftwnd_create_plan(_nd, _nxyz, dir, FFTW_ESTIMATE);
if (pnd == NULL)
throw AllocationError("FFTWServerPlan: failed to create plan (pnd) !");
}
}
FFTWServerPlan::~FFTWServerPlan()
......@@ -85,11 +130,16 @@ FFTWServerPlan::Recreate(int n)
if (p) {
fftw_destroy_plan(p);
p = fftw_create_plan(n, _dir, FFTW_ESTIMATE);
if (p == NULL)
throw AllocationError("FFTWServerPlan::Recreate failed to create plan (p) !");
}
else {
rfftw_destroy_plan(rp);
rp = rfftw_create_plan(n, _dir, FFTW_ESTIMATE);
if (rp == NULL)
throw AllocationError("FFTWServerPlan::Recreate failed to create plan (rp) !");
}
}
void
......@@ -116,10 +166,14 @@ FFTWServerPlan::Recreate(int nd, int * nxyz)
if (pnd) {
fftwnd_destroy_plan(pnd);
pnd = fftwnd_create_plan(_nd, _nxyz, _dir, FFTW_ESTIMATE);
if (pnd == NULL)
throw AllocationError("FFTWServerPlan::Recreate failed to create plan (pnd) !");
}
else {
rfftwnd_destroy_plan(rpnd);
rpnd = rfftwnd_create_plan(_nd, _nxyz, _dir, FFTW_ESTIMATE);
if (rpnd == NULL)
throw AllocationError("FFTWServerPlan::Recreate failed to create plan (rpnd) !");
}
}
......
......@@ -5,8 +5,10 @@
#include "fftservintf.h"
// Classe definissant l'interface pour les transformees de Fourier
// L'implementation par defaut est vide et lance une exception
// Classe implementant l'interface FFTServerInterface en
// utilisant FFTW
namespace SOPHYA {
class FFTWServerPlan;
......@@ -43,4 +45,6 @@ class FFTWServer : public FFTServerInterface {
FFTArrayChecker<r_8> ckR8;
};
} // Fin du namespace
#endif
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment