]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/TestZeroSuppress.C
Minor fixes
[u/mrichter/AliRoot.git] / FMD / scripts / TestZeroSuppress.C
1 TGraph*
2 makeInput(Int_t n, TArrayI& vals, TF1* f)
3 {
4   vals.Set(n);
5   TGraph* data = new TGraph(n);
6   data->SetName("input");
7   data->SetTitle("Input");
8   data->SetMarkerColor(kRed);
9   data->SetLineColor(kRed);
10   data->SetMarkerStyle(20);
11   
12   for (Int_t i = 0; i < data->GetN(); i++) { 
13     Double_t v = f->Eval(i);
14     data->SetPoint(i, i, v);
15     vals[i]   = Int_t(v);
16   }
17   return data;
18 }
19
20 TGraph* makeFlat(const char* name, Int_t n, Float_t val, TArrayF& vals)
21 {
22   vals.Set(n);
23   vals.Reset(val);
24   TGraph* g = new TGraph(n);
25   g->SetName(name);
26   g->SetTitle(name);
27   g->SetMarkerStyle(21);
28   
29   for (Int_t i = 0; i < g->GetN(); i++) 
30     g->SetPoint(i, i, val);
31   return g;
32 }
33   
34 void
35 runTest(TF1* f, AliFMDRawWriter* w, Float_t pv, Float_t nv)
36 {
37   static Int_t num = 1;
38   
39   Int_t    n         = Int_t(f->GetXmax());
40   UShort_t threshold = 1;
41   TArrayI  vals;
42   TArrayF  peds;
43   TArrayF  noise;
44   TArrayF  dummy;
45   TGraph*  input      = makeInput(n, vals, f);
46   TGraph*  gPeds      = makeFlat("peds",      n, pv,        peds);
47   TGraph*  gNoise     = makeFlat("noise",     n, nv,        noise);
48   TGraph*  gThreshold = makeFlat("threshold", n, threshold, dummy);
49   gThreshold->SetLineStyle(2);
50   gPeds->SetLineColor(kGreen);
51   gNoise->SetLineColor(kMagenta);
52
53   w->ZeroSuppress(vals.fArray, n, peds.fArray, noise.fArray, threshold);
54
55   TGraph* output = new TGraph(n);
56   output->SetName("output");
57   output->SetLineColor(kBlue);
58   output->SetMarkerColor(kBlue);
59   output->SetMarkerStyle(21);
60   for (Int_t i = 0; i < output->GetN(); i++) 
61     output->SetPoint(i, i, vals[i]);
62
63   TCanvas* c = new TCanvas(Form("c%02d", num), 
64                            Form("Zero suppression test %d", num));
65   c->SetFillColor(kWhite);
66   input->Draw("APL");
67   input->GetHistogram()->SetMinimum(0);
68   gPeds->Draw("L same");
69   gThreshold->Draw("L same");
70   output->Draw("PL same");
71   f->SetLineStyle(3);
72   f->Draw("same");
73   num++;
74 }
75
76   
77
78 void
79 TestZeroSuppress()
80 {
81   AliLog::SetModuleDebugLevel("FMD", 1);
82   AliFMDRawWriter* w = new AliFMDRawWriter(0);
83
84   TF1* simple = new TF1("simple", "[0] + [1] * (x >= [2] && x <= [3])", 0, 16);
85   simple->SetParameters(4, 10, 2.5, 6.5);
86
87   runTest(simple, w, 3.5, 0);
88
89   w.SetPedSubtract();
90   runTest(simple, w, 3.5, 0);
91
92   w.SetPedSubtract(kFALSE);
93   simple->SetParameters(4, 10, 3.5, 4.5);
94   runTest(simple, w, 3.5, 0);
95
96   TF1* two = new TF1("two", 
97                      "[0]+[1]*((x>=[2]&&x<=[3])||(x>=[4]&&x<=[5]))", 0, 16);
98   two->SetParameters(4, 10, 2.5, 4.5, 9.5, 12.5);
99   runTest(two, w, 3.5, 0);
100
101 }
102
103