FairRoot/PandaRoot
Public Member Functions | Private Member Functions | Private Attributes | List of all members
PndGemFindHitsAna Class Reference

analyze found GEM hits More...

#include <PndGemFindHitsAna.h>

Inheritance diagram for PndGemFindHitsAna:

Public Member Functions

 PndGemFindHitsAna ()
 
 PndGemFindHitsAna (Int_t iVerbose)
 
 PndGemFindHitsAna (TString taskName, Int_t iVerbose)
 
virtual ~PndGemFindHitsAna ()
 
virtual void Exec (Option_t *opt)
 
void SetVerbose (const Int_t &verbose)
 

Private Member Functions

void CreateHistos ()
 
void AnaHistos ()
 
virtual void SetParContainers ()
 
virtual void Finish ()
 
virtual InitStatus Init ()
 
virtual InitStatus ReInit ()
 
 ClassDef (PndGemFindHitsAna, 1)
 

Private Attributes

PndGemDigiParfDigiPar
 
TClonesArray * fGemHitArray
 
Int_t fNofEvents
 event counter More...
 
Double_t fGridSize
 
std::vector< Int_t > fStatBegHist
 
std::vector< Int_t > fGridHalfLen
 
TList * fHistoList
 
TClonesArray * fhFrontBackDiff
 

Detailed Description

analyze found GEM hits

Author
R. Karabowicz r.kar.nosp@m.abow.nosp@m.icz@g.nosp@m.si.d.nosp@m.e
Date
19.03.2009 Compare the found hits in all stations. Draw the ( position back - position front ). On each station two hits are registered, the difference in Z position of drift volumes are about 3 cm. Thus the difference in hits X,Y position should be rather well defined. From this it might be possible to align planes in each station.

Definition at line 32 of file PndGemFindHitsAna.h.

Constructor & Destructor Documentation

PndGemFindHitsAna::PndGemFindHitsAna ( )

Default constructor

Definition at line 34 of file PndGemFindHitsAna.cxx.

35  : FairTask("GEM Find Hits Ana", 1),
36  fDigiPar(NULL),
37  fGemHitArray(NULL),
38  fNofEvents(0),
39  fGridSize(10.),
40  // fStatBegHist(NULL),
41  // fGridHalfLen(NULL),
42  fHistoList(NULL),
43  fhFrontBackDiff(NULL)
44 {
45 }
PndGemDigiPar * fDigiPar
Int_t fNofEvents
event counter
TClonesArray * fhFrontBackDiff
TClonesArray * fGemHitArray
PndGemFindHitsAna::PndGemFindHitsAna ( Int_t  iVerbose)

Definition at line 49 of file PndGemFindHitsAna.cxx.

50  : FairTask("GEM Find Hits Ana", iVerbose) ,
51  fDigiPar(NULL),
52  fGemHitArray(NULL),
53  fNofEvents(0),
54  fGridSize(10.),
55  // fStatBegHist(NULL),
56  // fGridHalfLen(NULL),
57  fHistoList(NULL),
58  fhFrontBackDiff(NULL)
59 {
60 }
PndGemDigiPar * fDigiPar
Int_t fNofEvents
event counter
TClonesArray * fhFrontBackDiff
TClonesArray * fGemHitArray
Int_t iVerbose
PndGemFindHitsAna::PndGemFindHitsAna ( TString  taskName,
Int_t  iVerbose 
)

Definition at line 64 of file PndGemFindHitsAna.cxx.

65  : FairTask(taskName, iVerbose) ,
66  fDigiPar(NULL),
67  fGemHitArray(NULL),
68  fNofEvents(0),
69  fGridSize(10.),
70  // fStatBegHist(NULL),
71  // fGridHalfLen(NULL),
72  fHistoList(NULL),
73  fhFrontBackDiff(NULL)
74 {
75 }
PndGemDigiPar * fDigiPar
Int_t fNofEvents
event counter
TClonesArray * fhFrontBackDiff
TClonesArray * fGemHitArray
Int_t iVerbose
PndGemFindHitsAna::~PndGemFindHitsAna ( )
virtual

Destructor

Definition at line 79 of file PndGemFindHitsAna.cxx.

79 { }

Member Function Documentation

void PndGemFindHitsAna::AnaHistos ( )
private

Definition at line 273 of file PndGemFindHitsAna.cxx.

References Double_t, fDigiPar, fGridHalfLen, fGridSize, fhFrontBackDiff, Fill(), fStatBegHist, and PndGemDigiPar::GetNStations().

Referenced by Finish().

273  {
274 
275  Int_t nofStations = fDigiPar->GetNStations();
276 
277  Double_t quartX[4] = {1.,-1., 1.,-1.};
278  Double_t quartY[4] = {1., 1.,-1.,-1.};
279  for ( Int_t istat = 0 ; istat < nofStations ; istat++ ) {
280 
281  Int_t histCounter = fStatBegHist[istat];
282 
283  Int_t histogramToFill = fStatBegHist[istat] + 4*fGridHalfLen[istat]*fGridHalfLen[istat];
284 
285  for ( Int_t iquart = 0 ; iquart < 4 ; iquart++ ) {
286  for ( Int_t igry = 0 ; igry < fGridHalfLen[istat] ; igry++ ) {
287  for ( Int_t igrx = 0 ; igrx < fGridHalfLen[istat] ; igrx++ ) {
288  if ( ((TH2F*)(fhFrontBackDiff->At(histCounter)))->GetEntries() > 0 ) {
289  ((TH2F*)fhFrontBackDiff->At(histogramToFill+0))->Fill(quartX[iquart]*(fGridSize/2.+fGridSize*igrx),
290  quartY[iquart]*(fGridSize/2.+fGridSize*igry),
291  ((TH2F*)fhFrontBackDiff->At(histCounter))->GetMean(1));
292  ((TH2F*)fhFrontBackDiff->At(histogramToFill+1))->Fill(quartX[iquart]*(fGridSize/2.+fGridSize*igrx),
293  quartY[iquart]*(fGridSize/2.+fGridSize*igry),
294  ((TH2F*)fhFrontBackDiff->At(histCounter))->GetMean(2));
295  }
296  else {
297  ((TH2F*)fhFrontBackDiff->At(histogramToFill+0))->Fill(quartX[iquart]*(fGridSize/2.+fGridSize*igrx),
298  quartY[iquart]*(fGridSize/2.+fGridSize*igry),
299  -1111.);
300  ((TH2F*)fhFrontBackDiff->At(histogramToFill+1))->Fill(quartX[iquart]*(fGridSize/2.+fGridSize*igrx),
301  quartY[iquart]*(fGridSize/2.+fGridSize*igry),
302  -1111.);
303  }
304  histCounter++;
305  }
306  }
307  }
308  }
309 
310 }
PndGemDigiPar * fDigiPar
Double_t
TClonesArray * fhFrontBackDiff
std::vector< Int_t > fStatBegHist
Int_t GetNStations()
Definition: PndGemDigiPar.h:45
std::vector< Int_t > fGridHalfLen
h1 Fill(hit_theta-point_theta)
cout<<"the Event No is "<< i<< endl;{{if(hit_array->GetEntriesFast()!=mc_array->GetEntriesFast()) continue;PndSdsHit *hit=(PndSdsHit *) hit_array-> At(j)
Definition: anaLmdCluster.C:71
PndGemFindHitsAna::ClassDef ( PndGemFindHitsAna  ,
 
)
private
void PndGemFindHitsAna::CreateHistos ( )
private

Definition at line 208 of file PndGemFindHitsAna.cxx.

References fDigiPar, fGridHalfLen, fGridSize, fhFrontBackDiff, fHistoList, fStatBegHist, fVerbose, PndGemDigiPar::GetNStations(), PndGemSensor::GetOuterRadius(), PndGemDigiPar::GetSensor(), PndGemSensor::GetZ0(), and sensor.

Referenced by Init().

208  {
209  fHistoList = new TList();
210 
211  if ( fVerbose )
212  cout << "ought to create histogram for " << fDigiPar->GetNStations() << " stations." << endl;
213 
214  Int_t nofStations = fDigiPar->GetNStations();
215  Int_t nofHists = 0;
216 
217  fhFrontBackDiff = new TClonesArray("TH2F", 5000);
218 
219  for ( Int_t istat = 0 ; istat < nofStations ; istat++ ) {
220  //PndGemStation* station = (PndGemStation*)fDigiPar->GetStation(istat);
222 
223  if ( fVerbose )
224  cout << "sensor @ z = " << sensor->GetZ0() << " has an outer radius of " << sensor->GetOuterRadius() << endl;
225 
226  Int_t gridHalfLen = TMath::Ceil(sensor->GetOuterRadius()/fGridSize);
227  fStatBegHist.push_back(nofHists);
228  fGridHalfLen.push_back(gridHalfLen);
229 
230  if ( fVerbose )
231  cout << " --> the grid of " << fGridSize << "x" << fGridSize << "cm^2 requires " << 4*gridHalfLen*gridHalfLen << " histograms" << endl;
232 
233  for ( Int_t iquart = 0 ; iquart < 4 ; iquart++ ) {
234  for ( Int_t igry = 0 ; igry < gridHalfLen ; igry++ ) {
235  for ( Int_t igrx = 0 ; igrx < gridHalfLen ; igrx++ ) {
236  new ((*fhFrontBackDiff)[nofHists++]) TH2F(Form("fhFrontBackCorr_s%d_q%d_x%d_y%d",istat,iquart,igrx,igry),
237  Form("front back corr on station %d, x(%s%.1f:%s%.1f), y(%s%.1f:%s%.1f)",istat,
238  (iquart%2?"-":""),igrx*fGridSize,(iquart%2?"-":""),(igrx+1.)*fGridSize,
239  (iquart/2?"-":""),igry*fGridSize,(iquart/2?"-":""),(igry+1.)*fGridSize),
240  200,-1.,1.,
241  200,-1.,1.);
242  if ( fVerbose )
243  cout << nofHists-1 << ": histogram " << fhFrontBackDiff->At(nofHists-1)->GetName() << " / " << fhFrontBackDiff->At(nofHists-1)->GetTitle() << endl;
244  }
245  }
246  }
247 
248  new ((*fhFrontBackDiff)[nofHists++]) TH2F(Form("fhFrontBackDiffX_s%d",istat),
249  Form("front back X difference on station %d",istat),
250  2*gridHalfLen,-gridHalfLen*fGridSize,gridHalfLen*fGridSize,
251  2*gridHalfLen,-gridHalfLen*fGridSize,gridHalfLen*fGridSize);
252  new ((*fhFrontBackDiff)[nofHists++]) TH2F(Form("fhFrontBackDiffY_s%d",istat),
253  Form("front back Y difference on station %d",istat),
254  2*gridHalfLen,-gridHalfLen*fGridSize,gridHalfLen*fGridSize,
255  2*gridHalfLen,-gridHalfLen*fGridSize,gridHalfLen*fGridSize);
256  }
257 
258  if ( fVerbose ) {
259  cout << "HM, seems that fhFrontBackDiff is filled with histograms, and there are " << nofHists << " of them" << endl;
260 
261  cout << "fStatBegHist has " << fStatBegHist.size() << " entries:" << endl;
262  for ( size_t isbh = 0 ; isbh < fStatBegHist.size() ; isbh++ )
263  cout << " ---> " << fStatBegHist[isbh] << endl;
264  cout << "fGridHalfLen has " << fGridHalfLen.size() << " entries:" << endl;
265  for ( size_t ighl = 0 ; ighl < fGridHalfLen.size() ; ighl++ )
266  cout << " ---> " << fGridHalfLen[ighl] << endl;
267  }
268 
269 }
int fVerbose
Definition: poormantracks.C:24
PndGemDigiPar * fDigiPar
TGeoVolume * sensor
PndGemSensor * GetSensor(Int_t stationNr, Int_t sensorNr)
TClonesArray * fhFrontBackDiff
std::vector< Int_t > fStatBegHist
Int_t GetNStations()
Definition: PndGemDigiPar.h:45
std::vector< Int_t > fGridHalfLen
Double_t GetZ0() const
Definition: PndGemSensor.h:100
Double_t GetOuterRadius() const
Definition: PndGemSensor.h:103
void PndGemFindHitsAna::Exec ( Option_t *  opt)
virtual

Execution

Definition at line 153 of file PndGemFindHitsAna.cxx.

References CAMath::Abs(), fGemHitArray, fGridHalfLen, fGridSize, fhFrontBackDiff, Fill(), fNofEvents, fStatBegHist, fVerbose, and PndGemHit::GetStationNr().

153  {
154 
155  if ( fVerbose ) {
156  cout << "IN EVENT GOT " << fGemHitArray->GetEntries() << " hits." << endl;
157  cout << " - - - - - > " << fhFrontBackDiff->GetEntries() << " histograms" << endl;
158  }
159 
160  Int_t nofGemHits = fGemHitArray->GetEntries();
161 
162  PndGemHit* gemHit1 = NULL;
163  PndGemHit* gemHit2 = NULL;
164 
165  for ( Int_t igem1 = 0 ; igem1 < nofGemHits ; igem1++ ) {
166  gemHit1 = (PndGemHit*)fGemHitArray->At(igem1);
167  for ( Int_t igem2 = 0 ; igem2 < nofGemHits ; igem2++ ) {
168  gemHit2 = (PndGemHit*)fGemHitArray->At(igem2);
169 
170  if ( TMath::Abs( gemHit2->GetZ() - gemHit1->GetZ() - 2. ) < 2. ) {
171  Int_t quartId = 0;
172  if ( gemHit1->GetY() < 0. ) quartId += 2;
173  if ( gemHit1->GetX() < 0. ) quartId += 1;
174  Int_t histogramToFill =
175  fStatBegHist[gemHit1->GetStationNr()-1] +
176  fGridHalfLen[gemHit1->GetStationNr()-1]*fGridHalfLen[gemHit1->GetStationNr()-1]*quartId +
177  fGridHalfLen[gemHit1->GetStationNr()-1]*(Int_t(TMath::Abs(gemHit1->GetY())/fGridSize)) +
178  (Int_t(TMath::Abs(gemHit1->GetX())/fGridSize));
179 
180  ((TH2F*)fhFrontBackDiff->At(histogramToFill))->Fill(gemHit2->GetX()-gemHit1->GetX(),
181  gemHit2->GetY()-gemHit1->GetY());
182 
183  if ( fVerbose > 1 ) {
184  cout << "will compare " << gemHit1->GetX() << " " << gemHit1->GetY() << " " << gemHit1->GetZ() << endl
185  << " with " << gemHit2->GetX() << " " << gemHit2->GetY() << " " << gemHit2->GetZ() << endl;
186  cout << " IN STATION " << gemHit1->GetStationNr() << endl;
187  cout << " 1 : " << fStatBegHist[gemHit1->GetStationNr()-1] << endl;
188  cout << " 2 : " << fGridHalfLen[gemHit1->GetStationNr()-1]*fGridHalfLen[gemHit1->GetStationNr()-1]*quartId << endl;
189  cout << " 3 : " << fGridHalfLen[gemHit1->GetStationNr()-1]*(Int_t(gemHit1->GetY()/fGridSize)) << endl;
190  cout << " 4 : " << (Int_t(gemHit1->GetX()/fGridSize)) << endl;
191  cout << " histogram to fill " << histogramToFill << endl;
192  cout << " - - - - - > "
193  << fhFrontBackDiff->At(histogramToFill)->GetName() << " / \""
194  << fhFrontBackDiff->At(histogramToFill)->GetTitle() << "\"" << endl;
195  cout << "====================================================================================" << endl;
196  }
197  }
198  }
199  }
200 
201  fNofEvents++;
202 
203 }
int fVerbose
Definition: poormantracks.C:24
Int_t fNofEvents
event counter
static T Abs(const T &x)
Definition: PndCAMath.h:39
TClonesArray * fhFrontBackDiff
std::vector< Int_t > fStatBegHist
Int_t GetStationNr() const
Definition: PndGemHit.h:81
TClonesArray * fGemHitArray
std::vector< Int_t > fGridHalfLen
h1 Fill(hit_theta-point_theta)
void PndGemFindHitsAna::Finish ( )
privatevirtual

Finish

Definition at line 314 of file PndGemFindHitsAna.cxx.

References AnaHistos(), fhFrontBackDiff, fHistoList, fNofEvents, and next.

314  {
315  cout << "-------------------- PndGemFindHitsAna : Summary ------------------" << endl;
316  cout << " Events: " << setw(10) << fNofEvents << endl;
317  cout << "---------------------------------------------------------------------" << endl;
318 
319  AnaHistos();
320 
321  TFile* temp = gFile;
322  FairRootManager* ioman = FairRootManager::Instance();
323  gFile = ioman->GetOutFile();
324  gDirectory = (TDirectory*)gFile;
325 
326  gDirectory->mkdir("GemFindHitsAna");
327  gDirectory->cd("GemFindHitsAna");
328  TIter next(fHistoList);
329  for ( Int_t ihist = 0 ; ihist < fhFrontBackDiff->GetEntries() ; ihist++ ) {
330  fhFrontBackDiff->At(ihist)->Write();
331  }
332  // while ( TH1* histo = ((TH1*)next()) ) histo->Write();
333  gDirectory->cd("..");
334 
335  gFile = temp;
336 
337 }
Int_t fNofEvents
event counter
TClonesArray * fhFrontBackDiff
static int next[96]
Definition: ranlxd.cxx:374
InitStatus PndGemFindHitsAna::Init ( )
privatevirtual

Intialisation

Definition at line 82 of file PndGemFindHitsAna.cxx.

References CreateHistos(), fDigiPar, fGemHitArray, fGridSize, PndGemDigiPar::GetNStations(), and rtdb.

82  {
83 
84  // Get and check FairRootManager
85  FairRootManager* ioman = FairRootManager::Instance();
86  if( !ioman ) {
87  cout << "-E- "<< GetName() <<"::Init: "
88  << "RootManager not instantised!" << endl;
89  return kERROR;
90  }
91 
92  // Get the pointer to the singleton FairRunAna object
93  FairRunAna* ana = FairRunAna::Instance();
94  if(NULL == ana) {
95  cout << "-E- "<< GetName() <<"::Init :"
96  <<" no FairRunAna object!" << endl;
97  return kERROR;
98  }
99  // Get the pointer to run-time data base
100  FairRuntimeDb* rtdb = ana->GetRuntimeDb();
101  if(NULL == rtdb) {
102  cout << "-E- "<< GetName() <<"::Init :"
103  <<" no runtime database!" << endl;
104  return kERROR;
105  }
106 
107  // Get Gem hit Array
108  fGemHitArray = (TClonesArray*) ioman->GetObject("GEMHit");
109  if ( !fGemHitArray ) {
110  cout << "-E- " << GetName() << "::Init: No PndGemHit array!" << endl;
111  return kERROR;
112  }
113 
114  // Get GEM digitisation parameter container
115  fDigiPar = (PndGemDigiPar*)(rtdb->getContainer("PndGemDetectors"));
116 
117  cout << "-I- " << fName.Data() << "::Init(). There are " << fDigiPar->GetNStations() << " GEM stations." << endl;
118  cout << "-I- " << fName.Data() << "::Init(). Initialization succesfull." << endl;
119 
120  fGridSize = 5.;
121  CreateHistos();
122 
123  return kSUCCESS;
124 }
PndGemDigiPar * fDigiPar
Digitization Parameter Class for GEM part.
Definition: PndGemDigiPar.h:31
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
Int_t GetNStations()
Definition: PndGemDigiPar.h:45
TClonesArray * fGemHitArray
InitStatus PndGemFindHitsAna::ReInit ( )
privatevirtual

Reinitialisation

Definition at line 128 of file PndGemFindHitsAna.cxx.

128  {
129 
130  return kSUCCESS;
131 }
void PndGemFindHitsAna::SetParContainers ( )
privatevirtual

Get parameter containers

Definition at line 135 of file PndGemFindHitsAna.cxx.

References fDigiPar, fNofEvents, and run.

135  {
136 
137  // Get run and runtime database
138  FairRunAna* run = FairRunAna::Instance();
139  if ( ! run ) Fatal("SetParContainers", "No analysis run");
140 
141  FairRuntimeDb* db = run->GetRuntimeDb();
142  if ( ! db ) Fatal("SetParContainers", "No runtime database");
143 
144  // Get GEM digitisation parameter container
145  fDigiPar = (PndGemDigiPar*)(db->getContainer("PndGemDetectors"));
146 
147  fNofEvents = 0;
148 }
Int_t run
Definition: autocutx.C:47
PndGemDigiPar * fDigiPar
Int_t fNofEvents
event counter
Digitization Parameter Class for GEM part.
Definition: PndGemDigiPar.h:31
void PndGemFindHitsAna::SetVerbose ( const Int_t &  verbose)
inline

Public modifiers

Definition at line 50 of file PndGemFindHitsAna.h.

References fVerbose, and verbose.

50 { fVerbose = verbose; };
int fVerbose
Definition: poormantracks.C:24
#define verbose

Member Data Documentation

PndGemDigiPar* PndGemFindHitsAna::fDigiPar
private

Definition at line 50 of file PndGemFindHitsAna.h.

Referenced by AnaHistos(), CreateHistos(), Init(), and SetParContainers().

TClonesArray* PndGemFindHitsAna::fGemHitArray
private

Definition at line 57 of file PndGemFindHitsAna.h.

Referenced by Exec(), and Init().

std::vector<Int_t> PndGemFindHitsAna::fGridHalfLen
private

Definition at line 65 of file PndGemFindHitsAna.h.

Referenced by AnaHistos(), CreateHistos(), and Exec().

Double_t PndGemFindHitsAna::fGridSize
private

Definition at line 62 of file PndGemFindHitsAna.h.

Referenced by AnaHistos(), CreateHistos(), Exec(), and Init().

TClonesArray* PndGemFindHitsAna::fhFrontBackDiff
private

Definition at line 70 of file PndGemFindHitsAna.h.

Referenced by AnaHistos(), CreateHistos(), Exec(), and Finish().

TList* PndGemFindHitsAna::fHistoList
private

Definition at line 67 of file PndGemFindHitsAna.h.

Referenced by CreateHistos(), and Finish().

Int_t PndGemFindHitsAna::fNofEvents
private

event counter

Event counter

Definition at line 60 of file PndGemFindHitsAna.h.

Referenced by Exec(), Finish(), and SetParContainers().

std::vector<Int_t> PndGemFindHitsAna::fStatBegHist
private

Definition at line 64 of file PndGemFindHitsAna.h.

Referenced by AnaHistos(), CreateHistos(), and Exec().


The documentation for this class was generated from the following files: