]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCCalROC.cxx
Removing inheritances from TAttLine, TAttMarker and AliRndm in AliModule. The copy...
[u/mrichter/AliRoot.git] / TPC / AliTPCCalROC.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16
17 ///////////////////////////////////////////////////////////////////////////////
18 //                                                                           //
19 //  Calibration base class for a single ROC                                  //
20 //  Contains one float value per pad                                         //
21 //     mapping of the pads taken form AliTPCROC                              //
22 //                                                                           //
23 ///////////////////////////////////////////////////////////////////////////////
24
25 #include "AliTPCCalROC.h"
26 #include "TMath.h"
27 #include "TClass.h"
28 #include "TFile.h"
29 #include "TH2F.h"
30 ClassImp(AliTPCCalROC)
31
32
33 //_____________________________________________________________________________
34 AliTPCCalROC::AliTPCCalROC()
35              :TObject(),
36               fSector(0),
37               fNChannels(0),
38               fNRows(0),
39               fIndexes(0),
40               fData(0)
41 {
42   //
43   // Default constructor
44   //
45
46 }
47
48 //_____________________________________________________________________________
49 AliTPCCalROC::AliTPCCalROC(UInt_t  sector)
50              :TObject(),
51               fSector(0),
52               fNChannels(0),
53               fNRows(0),
54               fIndexes(0),
55               fData(0)
56 {
57   //
58   // Constructor that initializes a given sector
59   //
60   fSector = sector;
61   fNChannels    =  AliTPCROC::Instance()->GetNChannels(fSector);
62   fNRows        =  AliTPCROC::Instance()->GetNRows(fSector);
63   fIndexes      =  AliTPCROC::Instance()->GetRowIndexes(fSector);
64   fData = new Float_t[fNChannels];
65   for (UInt_t  idata = 0; idata< fNChannels; idata++) fData[idata] = 0.;
66 }
67
68 //_____________________________________________________________________________
69 AliTPCCalROC::AliTPCCalROC(const AliTPCCalROC &c)
70              :TObject(c),
71               fSector(0),
72               fNChannels(0),
73               fNRows(0),
74               fIndexes(0),
75               fData(0)
76 {
77   //
78   // AliTPCCalROC copy constructor
79   //
80   fSector = c.fSector;
81   fNChannels    =  AliTPCROC::Instance()->GetNChannels(fSector);
82   fNRows        =  AliTPCROC::Instance()->GetNRows(fSector);
83   fIndexes      =  AliTPCROC::Instance()->GetRowIndexes(fSector);
84   //
85   fData   = new Float_t[fNChannels];
86   for (UInt_t  idata = 0; idata< fNChannels; idata++) fData[idata] = c.fData[idata];
87 }
88 //____________________________________________________________________________
89 AliTPCCalROC & AliTPCCalROC::operator =(const AliTPCCalROC & param)
90 {
91   //
92   // assignment operator - dummy
93   //
94   fData=param.fData;
95   return (*this);
96 }
97
98
99 //_____________________________________________________________________________
100 AliTPCCalROC::~AliTPCCalROC()
101 {
102   //
103   // AliTPCCalROC destructor
104   //
105   if (fData) {
106     delete [] fData;
107     fData = 0;
108   }
109 }
110
111
112
113 void AliTPCCalROC::Streamer(TBuffer &R__b)
114 {
115    // Stream an object of class AliTPCCalROC.
116    if (R__b.IsReading()) {
117       AliTPCCalROC::Class()->ReadBuffer(R__b, this);
118       fIndexes =  AliTPCROC::Instance()->GetRowIndexes(fSector);
119    } else {
120       AliTPCCalROC::Class()->WriteBuffer(R__b,this);
121    }
122 }
123
124
125
126 void AliTPCCalROC::Draw(Option_t* option){
127   //
128   // create histogram with values and draw it
129   //
130   UInt_t maxPad = 0;
131   for (UInt_t irow=0; irow<fNRows; irow++){
132     if (GetNPads(irow)>maxPad) maxPad = GetNPads(irow);
133   }
134   char  name[1000];
135   sprintf(name,"%s ROC%d",GetTitle(),fSector);
136   TH2F * his = new TH2F(name,name,fNRows+10,-5, fNRows+5, maxPad+10, -(Int_t(maxPad/2))-5, maxPad/2+5);
137   for (UInt_t irow=0; irow<fNRows; irow++){
138     UInt_t npads = (Int_t)GetNPads(irow);
139     for (UInt_t ipad=0; ipad<=npads; ipad++){
140       his->Fill(irow+0.5,Int_t(ipad)-Int_t(npads/2)+0.5,GetValue(irow,ipad));
141     }
142   }
143   his->Draw(option);
144 }
145
146
147 void AliTPCCalROC::Test(){
148   //
149   // example function to show functionality and tes AliTPCCalROC
150   //
151   AliTPCCalROC  roc0(0);  
152   for (UInt_t irow = 0; irow <roc0.GetNrows(); irow++){
153     for (UInt_t ipad = 0; ipad <roc0.GetNPads(irow); ipad++){
154       Float_t value  = irow+ipad/1000.;
155       roc0.SetValue(irow,ipad,value);
156     }
157   }
158   //
159   AliTPCCalROC roc1(roc0);
160   for (UInt_t irow = 0; irow <roc1.GetNrows(); irow++){
161     for (UInt_t ipad = 0; ipad <roc1.GetNPads(irow); ipad++){
162       Float_t value  = irow+ipad/1000.;
163       if (roc1.GetValue(irow,ipad)!=value){
164         printf("Read/Write error\trow=%d\tpad=%d\n",irow,ipad);
165       }
166     }
167   }  
168   TFile f("calcTest.root","recreate");
169   roc0.Write("Roc0");
170   AliTPCCalROC * roc2 = (AliTPCCalROC*)f.Get("Roc0");
171   f.Close();
172   //
173   for (UInt_t irow = 0; irow <roc0.GetNrows(); irow++){
174     if (roc0.GetNPads(irow)!=roc2->GetNPads(irow))
175       printf("NPads - Read/Write error\trow=%d\n",irow);
176     for (UInt_t ipad = 0; ipad <roc1.GetNPads(irow); ipad++){
177       Float_t value  = irow+ipad/1000.;
178       if (roc2->GetValue(irow,ipad)!=value){
179         printf("Read/Write error\trow=%d\tpad=%d\n",irow,ipad);
180       }
181     }
182   }   
183 }
184