Project
Loading...
Searching...
No Matches
testO2TPCIDCFourierTransform.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
16
17#define BOOST_TEST_MODULE Test TPC O2TPCIDCFourierTransform class
18#define BOOST_TEST_MAIN
19#define BOOST_TEST_DYN_LINK
20#include <boost/test/unit_test.hpp>
22#include "TRandom.h"
23#include <numeric>
24
25namespace o2::tpc
26{
27
28static constexpr float ABSTOLERANCE = 0.01f; // absolute tolerance is taken at small values near 0
29static constexpr float TOLERANCE = 0.4f; // difference between original 1D-IDC and 1D-IDC from fourier transform -> inverse fourier transform
30
31o2::tpc::IDCOne get1DIDCs(const std::vector<unsigned int>& integrationIntervals)
32{
33 const unsigned int nIDCs = std::accumulate(integrationIntervals.begin(), integrationIntervals.end(), static_cast<unsigned int>(0));
34 o2::tpc::IDCOne idcsOut;
35 std::vector<float> idcs(nIDCs);
36 for (auto& val : idcs) {
37 val = gRandom->Gaus(0, 0.2);
38 }
39 idcsOut.mIDCOne = std::move(idcs);
40 return idcsOut;
41}
42
43std::vector<unsigned int> getIntegrationIntervalsPerTF(const unsigned int integrationIntervals, const unsigned int tfs)
44{
45 std::vector<unsigned int> intervals;
46 intervals.reserve(tfs);
47 for (unsigned int i = 0; i < tfs; ++i) {
48 const unsigned int additionalInterval = (i % 3) ? 1 : 0; // in each integration inerval are either 10 or 11 values when having 128 orbits per TF and 12 orbits integration length
49 intervals.emplace_back(integrationIntervals + additionalInterval);
50 }
51 return intervals;
52}
53
54// testing FT of aggregator
55BOOST_AUTO_TEST_CASE(IDCFourierTransformAggregator_test)
56{
57 const unsigned int integrationIntervals = 10; // number of integration intervals for first TF
58 const unsigned int tfs = 200; // number of aggregated TFs
59 const unsigned int rangeIDC = 200; // number of IDCs used to calculate the fourier coefficients
60 const unsigned int nFourierCoeff = rangeIDC + 2; // number of fourier coefficients which will be calculated/stored needs to be the maximum value to be able to perform IFT
62 gRandom->SetSeed(0);
63
64 for (int iType = 0; iType < 2; ++iType) {
65 const bool fft = iType == 0 ? false : true;
66 FtType::setFFT(fft);
67 FtType::setNThreads(2);
68
69 FtType idcFourierTransform{rangeIDC, nFourierCoeff};
70 const auto intervalsPerTF = getIntegrationIntervalsPerTF(integrationIntervals, tfs);
71 idcFourierTransform.setIDCs(get1DIDCs(intervalsPerTF), intervalsPerTF);
72 idcFourierTransform.setIDCs(get1DIDCs(intervalsPerTF), intervalsPerTF);
73 idcFourierTransform.calcFourierCoefficients(tfs);
74
75 const std::vector<unsigned int> offsetIndex = idcFourierTransform.getLastIntervals();
76 const auto idcOneExpanded = idcFourierTransform.getExpandedIDCOne();
77 const auto inverseFourier = idcFourierTransform.inverseFourierTransform();
78 for (unsigned int interval = 0; interval < idcFourierTransform.getNIntervals(); ++interval) {
79 for (unsigned int index = 0; index < rangeIDC; ++index) {
80 const float origIDCOne = idcOneExpanded[index + offsetIndex[interval]];
81 const float iFTIDCOne = inverseFourier[interval][index];
82 if (std::fabs(origIDCOne) < ABSTOLERANCE) {
83 BOOST_CHECK_SMALL(iFTIDCOne - origIDCOne, ABSTOLERANCE);
84 } else {
85 BOOST_CHECK_CLOSE(iFTIDCOne, origIDCOne, TOLERANCE);
86 }
87 }
88 }
89 }
90}
91
92// testing FT of EPN
93BOOST_AUTO_TEST_CASE(IDCFourierTransformEPN_test)
94{
95 const int nIter = 100; // number of iterations
96 const unsigned int integrationIntervals = 10; // number of integration intervals for first TF
97 const unsigned int rangeIDC = 200; // number of IDCs used to calculate the fourier coefficients
98 const unsigned int tfs = rangeIDC / integrationIntervals; // number of aggregated TFs (minimum number of 1D-IDCs obtained by get1DIDCs must be>rangeIDC)
99 const unsigned int nFourierCoeff = rangeIDC + 2; // number of fourier coefficients which will be calculated/stored needs to be the maximum value to be able to perform IFT
101 gRandom->SetSeed(0);
102
103 for (int iter = 0; iter < nIter; ++iter) {
104 for (int iType = 0; iType < 2; ++iType) {
105 const bool fft = iType == 0 ? false : true;
106 FtType::setFFT(fft);
107 FtType idcFourierTransform{rangeIDC, nFourierCoeff};
108 const auto intervalsPerTF = getIntegrationIntervalsPerTF(integrationIntervals, tfs);
109 idcFourierTransform.setIDCs(get1DIDCs(intervalsPerTF));
110 idcFourierTransform.calcFourierCoefficients();
111 const auto offsetIndex = idcFourierTransform.getLastIntervals();
112 const auto idcOneExpanded = idcFourierTransform.getExpandedIDCOne();
113 const auto inverseFourier = idcFourierTransform.inverseFourierTransform();
114 for (unsigned int interval = 0; interval < idcFourierTransform.getNIntervals(); ++interval) {
115 for (unsigned int index = 0; index < rangeIDC; ++index) {
116 const float origIDCOne = idcOneExpanded[index + offsetIndex[interval]];
117 const float iFTIDCOne = inverseFourier[interval][index];
118 if (std::fabs(origIDCOne) < ABSTOLERANCE) {
119 BOOST_CHECK_SMALL(iFTIDCOne - origIDCOne, ABSTOLERANCE);
120 } else {
121 BOOST_CHECK_CLOSE(iFTIDCOne, origIDCOne, TOLERANCE);
122 }
123 }
124 }
125 }
126 }
127}
128
129} // namespace o2::tpc
int32_t i
class for calculating the fourier coefficients from 1D-IDCs
GLuint index
Definition glcorearb.h:781
GLuint GLfloat * val
Definition glcorearb.h:1582
Global TPC definitions and constants.
Definition SimTraits.h:167
BOOST_AUTO_TEST_CASE(ClusterHardware_test1)
std::vector< unsigned int > getIntegrationIntervalsPerTF(const unsigned int integrationIntervals, const unsigned int tfs)
o2::tpc::IDCOne get1DIDCs(const std::vector< unsigned int > &integrationIntervals)
std::vector< float > mIDCOne
I_1(t) = <I(r,\phi,t) / I_0(r,\phi)>_{r,\phi}.