Skip to content
Snippets Groups Projects
Commit b7c90c9c authored by Mark OLESEN's avatar Mark OLESEN
Browse files

COMP: label-size 64 compilation of fft (issue #813)

parent ba1ad60a
No related branches found
No related tags found
No related merge requests found
...@@ -3,7 +3,7 @@ ...@@ -3,7 +3,7 @@
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | \\ / O peration |
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation | Copyright (C) 2016 OpenCFD Ltd. \\/ M anipulation | Copyright (C) 2016-2018 OpenCFD Ltd.
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
This file is part of OpenFOAM. This file is part of OpenFOAM.
...@@ -24,7 +24,6 @@ License ...@@ -24,7 +24,6 @@ License
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
#include "fft.H" #include "fft.H"
#include <fftw3.h> #include <fftw3.h>
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
...@@ -33,7 +32,7 @@ void Foam::fft::fftRenumberRecurse ...@@ -33,7 +32,7 @@ void Foam::fft::fftRenumberRecurse
( (
List<complex>& data, List<complex>& data,
List<complex>& renumData, List<complex>& renumData,
const labelList& nn, const UList<int>& nn,
label nnprod, label nnprod,
label ii, label ii,
label l1, label l1,
...@@ -89,7 +88,7 @@ void Foam::fft::fftRenumberRecurse ...@@ -89,7 +88,7 @@ void Foam::fft::fftRenumberRecurse
} }
void Foam::fft::fftRenumber(List<complex>& data, const labelList& nn) void Foam::fft::fftRenumber(List<complex>& data, const UList<int>& nn)
{ {
List<complex> renumData(data); List<complex> renumData(data);
...@@ -181,6 +180,8 @@ void Foam::fft::transform ...@@ -181,6 +180,8 @@ void Foam::fft::transform
) )
{ {
// Copy field into fftw containers // Copy field into fftw containers
// NB: this is not fully compliant, since array sizing is non-constexpr.
// However, cannot instantiate List with fftw_complex
const label N = field.size(); const label N = field.size();
fftw_complex in[N], out[N]; fftw_complex in[N], out[N];
...@@ -236,13 +237,13 @@ Foam::tmp<Foam::complexField> Foam::fft::forwardTransform ...@@ -236,13 +237,13 @@ Foam::tmp<Foam::complexField> Foam::fft::forwardTransform
const UList<int>& nn const UList<int>& nn
) )
{ {
tmp<complexField> tfftField(new complexField(tfield)); tmp<complexField> tresult(new complexField(tfield));
transform(tfftField.ref(), nn, FORWARD_TRANSFORM); transform(tresult.ref(), nn, FORWARD_TRANSFORM);
tfield.clear(); tfield.clear();
return tfftField; return tresult;
} }
...@@ -252,13 +253,13 @@ Foam::tmp<Foam::complexField> Foam::fft::reverseTransform ...@@ -252,13 +253,13 @@ Foam::tmp<Foam::complexField> Foam::fft::reverseTransform
const UList<int>& nn const UList<int>& nn
) )
{ {
tmp<complexField> tifftField(new complexField(tfield)); tmp<complexField> tresult(new complexField(tfield));
transform(tifftField.ref(), nn, REVERSE_TRANSFORM); transform(tresult.ref(), nn, REVERSE_TRANSFORM);
tfield.clear(); tfield.clear();
return tifftField; return tresult;
} }
...@@ -268,17 +269,11 @@ Foam::tmp<Foam::complexVectorField> Foam::fft::forwardTransform ...@@ -268,17 +269,11 @@ Foam::tmp<Foam::complexVectorField> Foam::fft::forwardTransform
const UList<int>& nn const UList<int>& nn
) )
{ {
tmp<complexVectorField> tfftVectorField auto tresult = tmp<complexVectorField>::New(tfield().size());
(
new complexVectorField
(
tfield().size()
)
);
for (direction cmpt=0; cmpt<vector::nComponents; cmpt++) for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
{ {
tfftVectorField.ref().replace tresult.ref().replace
( (
cmpt, cmpt,
forwardTransform(tfield().component(cmpt), nn) forwardTransform(tfield().component(cmpt), nn)
...@@ -287,7 +282,7 @@ Foam::tmp<Foam::complexVectorField> Foam::fft::forwardTransform ...@@ -287,7 +282,7 @@ Foam::tmp<Foam::complexVectorField> Foam::fft::forwardTransform
tfield.clear(); tfield.clear();
return tfftVectorField; return tresult;
} }
...@@ -297,17 +292,11 @@ Foam::tmp<Foam::complexVectorField> Foam::fft::reverseTransform ...@@ -297,17 +292,11 @@ Foam::tmp<Foam::complexVectorField> Foam::fft::reverseTransform
const UList<int>& nn const UList<int>& nn
) )
{ {
tmp<complexVectorField> tifftVectorField auto tresult = tmp<complexVectorField>::New(tfield().size());
(
new complexVectorField
(
tfield().size()
)
);
for (direction cmpt=0; cmpt<vector::nComponents; cmpt++) for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
{ {
tifftVectorField.ref().replace tresult.ref().replace
( (
cmpt, cmpt,
reverseTransform(tfield().component(cmpt), nn) reverseTransform(tfield().component(cmpt), nn)
...@@ -316,7 +305,7 @@ Foam::tmp<Foam::complexVectorField> Foam::fft::reverseTransform ...@@ -316,7 +305,7 @@ Foam::tmp<Foam::complexVectorField> Foam::fft::reverseTransform
tfield.clear(); tfield.clear();
return tifftVectorField; return tresult;
} }
......
...@@ -3,7 +3,7 @@ ...@@ -3,7 +3,7 @@
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | \\ / O peration |
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
\\/ M anipulation | Copyright (C) 2016 OpenCFD Ltd. \\/ M anipulation | Copyright (C) 2016-2018 OpenCFD Ltd.
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
This file is part of OpenFOAM. This file is part of OpenFOAM.
...@@ -57,15 +57,15 @@ public: ...@@ -57,15 +57,15 @@ public:
enum transformDirection enum transformDirection
{ {
FORWARD_TRANSFORM = -1, FORWARD_TRANSFORM = -1, //!< The sign -1 = FFTW_FORWARD
REVERSE_TRANSFORM = 1 REVERSE_TRANSFORM = 1 //!< The sign +1 = FFTW_BACKWARD
}; };
static void fftRenumberRecurse static void fftRenumberRecurse
( (
List<complex>& data, List<complex>& data,
List<complex>& renumData, List<complex>& renumData,
const labelList& nn, const UList<int>& nn,
label nnprod, label nnprod,
label ii, label ii,
label l1, label l1,
...@@ -74,7 +74,7 @@ public: ...@@ -74,7 +74,7 @@ public:
//- fftRenumber: fold the n-d data array to get the fft components in //- fftRenumber: fold the n-d data array to get the fft components in
//- the right places. //- the right places.
static void fftRenumber(List<complex>& data, const labelList& nn); static void fftRenumber(List<complex>& data, const UList<int>& nn);
//- Transform real-value data //- Transform real-value data
// - uses the fftw real to half-complex method // - uses the fftw real to half-complex method
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment