]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FORWARD/analysis/AliFMDAnaParameters.cxx
Upgrade following the 900 GeV data. A lot of extra things added - histograms etc...
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis / AliFMDAnaParameters.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 //The design of this class is based on the AliFMDParameters class. Its purpose
17 //is to hold parameters for the analysis such as background correction and 
18 //fit functions.
19 //
20 //Author: Hans Hjersing Dalsgaard, NBI, hans.dalsgaard@cern.ch
21 //
22
23 #include "AliFMDDebug.h"                   // ALILOG_H
24 #include "AliFMDAnaParameters.h"           // ALIFMDPARAMETERS_H
25 //#include <AliCDBManager.h>         // ALICDBMANAGER_H
26 //#include <AliCDBEntry.h>           // ALICDBMANAGER_H
27 //#include "AliFMDRing.h"
28 #include <AliLog.h>
29 #include <Riostream.h>
30 #include <sstream>
31 #include <TSystem.h>
32 #include <TH2D.h>
33 #include <TF1.h>
34 #include <TMath.h>
35 #include "AliTriggerAnalysis.h"
36 #include "AliPhysicsSelection.h"
37 #include "AliBackgroundSelection.h"
38 #include "AliMultiplicity.h"
39 #include "AliESDEvent.h"
40 #include "AliESDVertex.h"
41
42 //====================================================================
43 ClassImp(AliFMDAnaParameters)
44 #if 0
45   ; // This is here to keep Emacs for indenting the next line
46 #endif
47
48 //const char* AliFMDAnaParameters::fgkBackgroundCorrection  = "FMD/Correction/Background";
49 //const char* AliFMDAnaParameters::fgkEnergyDists = "FMD/Correction/EnergyDistribution";
50 const char* AliFMDAnaParameters::fgkBackgroundID         = "background";
51 const char* AliFMDAnaParameters::fgkEnergyDistributionID = "energydistributions";
52 const char* AliFMDAnaParameters::fgkEventSelectionEffID  = "eventselectionefficiency";
53 const char* AliFMDAnaParameters::fgkSharingEffID         = "sharingefficiency";
54 //____________________________________________________________________
55 AliFMDAnaParameters* AliFMDAnaParameters::fgInstance = 0;
56
57 //____________________________________________________________________
58
59 AliFMDAnaParameters* 
60 AliFMDAnaParameters::Instance() 
61 {
62   // Get static instance 
63   if (!fgInstance) fgInstance = new AliFMDAnaParameters;
64   return fgInstance;
65 }
66
67 //____________________________________________________________________
68 AliFMDAnaParameters::AliFMDAnaParameters() :
69   fIsInit(kFALSE),
70   fBackground(0),
71   fEnergyDistribution(0),
72   fEventSelectionEfficiency(0),
73   fSharingEfficiency(0),
74   fCorner1(4.2231, 26.6638),
75   fCorner2(1.8357, 27.9500),
76   fEnergyPath("$ALICE_ROOT/PWG2/FORWARD/corrections/EnergyDistribution"),
77   fBackgroundPath("$ALICE_ROOT/PWG2/FORWARD/corrections/Background"),
78   fEventSelectionEffPath("$ALICE_ROOT/PWG2/FORWARD/corrections/EventSelectionEfficiency"),
79   fSharingEffPath("$ALICE_ROOT/PWG2/FORWARD/corrections/SharingEfficiency"),
80   fProcessPrimary(kFALSE),
81   fProcessHits(kFALSE),
82   fTrigger(kMB1),
83   fEnergy(k10000),
84   fMagField(k5G),
85   fSpecies(kPP),
86   fPhysicsSelection(0),
87   fRealData(kFALSE),
88   fSPDlowLimit(0),
89   fSPDhighLimit(999999999),
90   fCentralSelection(kFALSE)
91 {
92   fPhysicsSelection = new AliPhysicsSelection;
93   AliBackgroundSelection* backgroundSelection = new AliBackgroundSelection("bg","bg");
94   
95   fPhysicsSelection->AddBackgroundIdentification(backgroundSelection);
96   //fPhysicsSelection->Initialize(104792);
97   // Do not use this - it is only for IO 
98   fgInstance = this;
99   // Default constructor 
100 }
101 //____________________________________________________________________
102 char* AliFMDAnaParameters::GetPath(const char* species) {
103   
104   char* path = "";
105   
106   if(species == fgkBackgroundID)
107     path = Form("%s/%s_%d_%d_%d_%d_%d_%d.root",
108                 fBackgroundPath.Data(),
109                 fgkBackgroundID,
110                 fEnergy,
111                 fTrigger,
112                 fMagField,
113                 fSpecies,
114                 0,
115                 0);
116   if(species == fgkEnergyDistributionID)
117     path = Form("%s/%s_%d_%d_%d_%d_%d_%d.root",
118                 fEnergyPath.Data(),
119                 fgkEnergyDistributionID,
120                 fEnergy,
121                 fTrigger,
122                 fMagField,
123                 fSpecies,
124                 0,
125                 0);
126   if(species == fgkEventSelectionEffID)
127     path = Form("%s/%s_%d_%d_%d_%d_%d_%d.root",
128                 fEventSelectionEffPath.Data(),
129                 fgkEventSelectionEffID,
130                 fEnergy,
131                 fTrigger,
132                 fMagField,
133                 fSpecies,
134                 0,
135                 0);
136   if(species == fgkSharingEffID)
137     path = Form("%s/%s_%d_%d_%d_%d_%d_%d.root",
138                 fSharingEffPath.Data(),
139                 fgkSharingEffID,
140                 fEnergy,
141                 fTrigger,
142                 fMagField,
143                 fSpecies,
144                 0,
145                 0);
146
147   return path;
148 }
149 //____________________________________________________________________
150 void AliFMDAnaParameters::Init(Bool_t forceReInit, UInt_t what)
151 {
152   // Initialize the parameters manager.  We need to get stuff from the
153   // CDB here. 
154   if (forceReInit) fIsInit = kFALSE;
155   if (fIsInit) return;
156   if (what & kBackgroundCorrection)       InitBackground();
157   if (what & kEnergyDistributions)        InitEnergyDists();
158   if (what & kEventSelectionEfficiency)   InitEventSelectionEff();
159   if (what & kSharingEfficiency)          InitSharingEff();
160   
161
162   
163   fIsInit = kTRUE;
164 }
165 //____________________________________________________________________
166
167 void AliFMDAnaParameters::InitBackground() {
168   
169   //AliCDBEntry*   background = GetEntry(fgkBackgroundCorrection);
170   
171   TFile* fin = TFile::Open(GetPath(fgkBackgroundID));
172   
173   if (!fin) return;
174   
175   fBackground = dynamic_cast<AliFMDAnaCalibBackgroundCorrection*>(fin->Get(fgkBackgroundID));
176   if (!fBackground) AliFatal("Invalid background object from CDB");
177   
178 }
179
180 //____________________________________________________________________
181
182 void AliFMDAnaParameters::InitEnergyDists() {
183   
184   TFile* fin = TFile::Open(GetPath(fgkEnergyDistributionID));
185   //AliCDBEntry*   edist = GetEntry(fgkEnergyDists);
186   if (!fin) return;
187   
188   fEnergyDistribution = dynamic_cast<AliFMDAnaCalibEnergyDistribution*>(fin->Get(fgkEnergyDistributionID));
189   
190   if (!fEnergyDistribution) AliFatal("Invalid background object from CDB");
191   
192 }
193
194 //____________________________________________________________________
195
196 void AliFMDAnaParameters::InitEventSelectionEff() {
197   
198   //AliCDBEntry*   background = GetEntry(fgkBackgroundCorrection);
199   TFile* fin = TFile::Open(GetPath(fgkEventSelectionEffID));
200                             
201   if (!fin) return;
202   
203   fEventSelectionEfficiency = dynamic_cast<AliFMDAnaCalibEventSelectionEfficiency*>(fin->Get(fgkEventSelectionEffID));
204   if (!fEventSelectionEfficiency) AliFatal("Invalid background object from CDB");
205   
206 }
207 //____________________________________________________________________
208
209 void AliFMDAnaParameters::PrintStatus() const
210 {
211   
212   TString energystring;
213   switch(fEnergy) {
214   case k900:
215     energystring.Form("900 GeV");   break;
216   case k7000:
217     energystring.Form("7000 GeV");  break;
218   case k10000:
219     energystring.Form("10000 GeV"); break;
220   case k14000:
221     energystring.Form("14000 GeV"); break;
222   default:
223     energystring.Form("invalid energy"); break;
224   }
225   TString triggerstring;
226   switch(fTrigger) {
227   case kMB1:
228     triggerstring.Form("Minimum bias 1");   break;
229   case kMB2:
230     triggerstring.Form("Minimum bias 2");   break;
231   case kSPDFASTOR:
232     triggerstring.Form("SPD FAST OR");   break;
233   case kNOCTP:
234     triggerstring.Form("NO TRIGGER TEST");   break;
235   default:
236     energystring.Form("invalid trigger"); break;
237   }
238   TString magstring;
239   switch(fMagField) {
240   case k5G:
241     magstring.Form("5 kGaus");   break;
242   case k0G:
243     magstring.Form("0 kGaus");   break;
244   default:
245     magstring.Form("invalid mag field"); break;
246   }
247   TString collsystemstring;
248   switch(fSpecies) {
249   case kPP:
250     collsystemstring.Form("p+p");   break;
251   case kPbPb:
252     collsystemstring.Form("Pb+Pb");   break;
253   default:
254     collsystemstring.Form("invalid collision system");   break;
255   }
256   
257
258   std::cout<<"Energy      = "<<energystring.Data()<<std::endl;
259   std::cout<<"Trigger     = "<<triggerstring.Data()<<std::endl;
260   std::cout<<"Mag Field   = "<<magstring.Data()<<std::endl;
261   std::cout<<"Coll System = "<<collsystemstring.Data()<<std::endl;
262   
263   
264   
265 }
266
267 //____________________________________________________________________
268
269 void AliFMDAnaParameters::InitSharingEff() {
270   
271   //AliCDBEntry*   background = GetEntry(fgkBackgroundCorrection);
272   TFile* fin = TFile::Open(GetPath(fgkSharingEffID));
273                             
274   if (!fin) return;
275   
276   fSharingEfficiency = dynamic_cast<AliFMDAnaCalibSharingEfficiency*>(fin->Get(fgkSharingEffID));
277   if (!fSharingEfficiency) AliFatal("Invalid background object from CDB");
278   
279 }
280 //____________________________________________________________________
281 Float_t AliFMDAnaParameters::GetVtxCutZ() {
282   
283   if(!fIsInit) {
284     AliWarning("Not initialized yet. Call Init() to remedy");
285     return -1;
286   }
287   
288   return fBackground->GetVtxCutZ();
289 }
290
291 //____________________________________________________________________
292 Int_t AliFMDAnaParameters::GetNvtxBins() {
293   
294   if(!fIsInit) {
295     AliWarning("Not initialized yet. Call Init() to remedy");
296     return -1;
297   }
298   
299   return fBackground->GetNvtxBins();
300 }
301 //____________________________________________________________________
302 TH1F* AliFMDAnaParameters::GetEnergyDistribution(Int_t det, Char_t ring, Float_t eta) {
303   
304   return fEnergyDistribution->GetEnergyDistribution(det, ring, eta);
305 }
306 //____________________________________________________________________
307 TH1F* AliFMDAnaParameters::GetEmptyEnergyDistribution(Int_t det, Char_t ring) {
308   
309   return fEnergyDistribution->GetEmptyEnergyDistribution(det, ring);
310 }
311 //____________________________________________________________________
312 TH1F* AliFMDAnaParameters::GetRingEnergyDistribution(Int_t det, Char_t ring) {
313   
314   return fEnergyDistribution->GetRingEnergyDistribution(det, ring);
315 }
316 //____________________________________________________________________
317 Float_t AliFMDAnaParameters::GetSigma(Int_t det, Char_t ring, Float_t eta) {
318   
319   if(!fIsInit) {
320     AliWarning("Not initialized yet. Call Init() to remedy");
321     return 0;
322   }
323   
324   TH1F* hEnergyDist       = GetEnergyDistribution(det,ring, eta);
325   TF1*  fitFunc           = hEnergyDist->GetFunction("FMDfitFunc");
326   if(!fitFunc) {
327     //AliWarning(Form("No function for FMD%d%c, eta %f",det,ring,eta));
328     return 1024;
329   }
330   Float_t sigma           = fitFunc->GetParameter(2);
331   return sigma;
332 }
333
334
335 //____________________________________________________________________
336 Float_t AliFMDAnaParameters::GetMPV(Int_t det, Char_t ring, Float_t eta) {
337   
338   if(!fIsInit) {
339     AliWarning("Not initialized yet. Call Init() to remedy");
340     return 0;
341   }
342   //AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
343   TH1F* hEnergyDist     = GetEnergyDistribution(det,ring,eta);
344   TF1*  fitFunc         = hEnergyDist->GetFunction("FMDfitFunc");
345   if(!fitFunc) {
346     AliWarning(Form("No function for FMD%d%c, eta %f",det,ring,eta));
347     std::cout<<hEnergyDist->GetEntries()<<"    "<<GetEtaBin(eta)<<std::endl;
348     return 1024;
349   }
350     
351   Float_t mpv           = fitFunc->GetParameter(1);
352   return mpv;
353 }
354 //____________________________________________________________________
355 Float_t AliFMDAnaParameters::GetConstant(Int_t det, Char_t ring, Float_t eta) {
356   
357   if(!fIsInit) {
358     AliWarning("Not initialized yet. Call Init() to remedy");
359     return 0;
360   }
361   
362   TH1F* hEnergyDist     = GetEnergyDistribution(det,ring,eta);
363   TF1*  fitFunc         = hEnergyDist->GetFunction("FMDfitFunc");
364   if(!fitFunc) {
365     AliWarning(Form("No function for FMD%d%c, eta %f",det,ring,eta));
366     return 0;
367   }
368     
369   Float_t mpv           = fitFunc->GetParameter(0);
370   return mpv;
371 }
372 //____________________________________________________________________
373 Float_t AliFMDAnaParameters::Get2MIPWeight(Int_t det, Char_t ring, Float_t eta) {
374   
375   if(!fIsInit) {
376     AliWarning("Not initialized yet. Call Init() to remedy");
377     return 0;
378   }
379   
380   TH1F* hEnergyDist     = GetEnergyDistribution(det,ring,eta);
381   TF1*  fitFunc         = hEnergyDist->GetFunction("FMDfitFunc");
382   if(!fitFunc) return 0;
383   Float_t twoMIPweight    = fitFunc->GetParameter(3);
384   
385   
386   
387   if(twoMIPweight < 1e-05)
388     twoMIPweight = 0;
389   
390   return twoMIPweight;
391 }
392 //____________________________________________________________________
393 Float_t AliFMDAnaParameters::Get3MIPWeight(Int_t det, Char_t ring, Float_t eta) {
394   
395   if(!fIsInit) {
396     AliWarning("Not initialized yet. Call Init() to remedy");
397     return 0;
398   }
399   
400   TH1F* hEnergyDist     = GetEnergyDistribution(det,ring,eta);
401   TF1*  fitFunc         = hEnergyDist->GetFunction("FMDfitFunc");
402   if(!fitFunc) return 0;
403   Float_t threeMIPweight    = fitFunc->GetParameter(4);
404   
405   if(threeMIPweight < 1e-05)
406     threeMIPweight = 0;
407   
408   Float_t twoMIPweight    = fitFunc->GetParameter(3);
409   
410   if(twoMIPweight < 1e-05)
411     threeMIPweight = 0;
412     
413   return threeMIPweight;
414 }
415 //____________________________________________________________________
416 Int_t AliFMDAnaParameters::GetNetaBins() {
417   return GetBackgroundCorrection(1,'I',0)->GetNbinsX();
418   
419 }
420 //____________________________________________________________________
421 Float_t AliFMDAnaParameters::GetEtaMin() {
422
423   return GetBackgroundCorrection(1,'I',0)->GetXaxis()->GetXmin();
424
425 //____________________________________________________________________
426 Float_t AliFMDAnaParameters::GetEtaMax() {
427
428 return GetBackgroundCorrection(1,'I',0)->GetXaxis()->GetXmax();
429
430 }
431 //____________________________________________________________________
432 Int_t AliFMDAnaParameters::GetEtaBin(Float_t eta) {
433   TAxis testaxis(GetNetaBins(),GetEtaMin(),GetEtaMax());
434   Int_t binnumber = testaxis.FindBin(eta) ;
435   
436   return binnumber;
437
438 }
439 //____________________________________________________________________
440
441 TH2F* AliFMDAnaParameters::GetBackgroundCorrection(Int_t det, 
442                                                    Char_t ring, 
443                                                    Int_t vtxbin) {
444   
445   if(!fIsInit) {
446     AliWarning("Not initialized yet. Call Init() to remedy");
447     return 0;
448   }
449   
450   
451   
452   if(vtxbin > fBackground->GetNvtxBins()) {
453     AliWarning(Form("No background object for vertex bin %d", vtxbin));
454     return 0;
455   } 
456   
457   return fBackground->GetBgCorrection(det,ring,vtxbin);
458 }
459 //____________________________________________________________________
460
461 TH1F* AliFMDAnaParameters::GetDoubleHitCorrection(Int_t det, 
462                                                   Char_t ring) {
463   
464   if(!fIsInit) {
465     AliWarning("Not initialized yet. Call Init() to remedy");
466     return 0;
467   }
468   
469   return fBackground->GetDoubleHitCorrection(det,ring);
470 }
471 //_____________________________________________________________________
472 Float_t AliFMDAnaParameters::GetEventSelectionEfficiency(Int_t vtxbin) {
473   if(!fIsInit) {
474     AliWarning("Not initialized yet. Call Init() to remedy");
475     return 0;
476   }
477   return fEventSelectionEfficiency->GetCorrection(vtxbin);
478
479 }
480 //_____________________________________________________________________
481 TH1F* AliFMDAnaParameters::GetSharingEfficiency(Int_t det, Char_t ring, Int_t vtxbin) {
482   if(!fIsInit) {
483     AliWarning("Not initialized yet. Call Init() to remedy");
484     return 0;
485   }
486   
487   return fSharingEfficiency->GetSharingEff(det,ring,vtxbin);
488
489 }
490 //_____________________________________________________________________
491 TH1F* AliFMDAnaParameters::GetSharingEfficiencyTrVtx(Int_t det, Char_t ring, Int_t vtxbin) {
492   if(!fIsInit) {
493     AliWarning("Not initialized yet. Call Init() to remedy");
494     return 0;
495   }
496   
497   return fSharingEfficiency->GetSharingEffTrVtx(det,ring,vtxbin);
498
499 }
500 //_____________________________________________________________________
501 Float_t AliFMDAnaParameters::GetMaxR(Char_t ring) const{
502   Float_t radius = 0;
503   if(ring == 'I')
504     radius = 17.2;
505   else if(ring == 'O')
506     radius = 28.0;
507   else
508     AliWarning("Unknown ring - must be I or O!");
509   
510   return radius;
511 }
512 //_____________________________________________________________________
513 Float_t AliFMDAnaParameters::GetMinR(Char_t ring) const{
514   Float_t radius = 0;
515   if(ring == 'I')
516     radius = 4.5213;
517   else if(ring == 'O')
518     radius = 15.4;
519   else
520     AliWarning("Unknown ring - must be I or O!");
521   
522   return radius;
523
524 }
525 //_____________________________________________________________________
526 void AliFMDAnaParameters::SetCorners(Char_t ring) {
527   
528   if(ring == 'I') {
529     fCorner1.Set(4.9895, 15.3560);
530     fCorner2.Set(1.8007, 17.2000);
531   }
532   else {
533     fCorner1.Set(4.2231, 26.6638);
534     fCorner2.Set(1.8357, 27.9500);
535   }
536   
537 }
538 //_____________________________________________________________________
539 Float_t AliFMDAnaParameters::GetPhiFromSector(UShort_t det, Char_t ring, UShort_t sec) const
540 {
541   Int_t nsec = (ring == 'I' ? 20 : 40);
542   Float_t basephi = 0;
543   if(det == 1) 
544     basephi = 1.72787594; 
545   if(det == 2 && ring == 'I')
546     basephi = 0.15707963;
547   if(det == 2 && ring == 'O')
548     basephi = 0.078539818;
549   if(det == 3 && ring == 'I')
550     basephi = 2.984513044;
551   if(det == 3 && ring == 'O')
552     basephi = 3.06305289;
553   
554   Float_t step = 2*TMath::Pi() / nsec;
555   Float_t phi = 0;
556   if(det == 3)
557     phi = basephi - sec*step;
558   else
559     phi = basephi + sec*step;
560   
561   if(phi < 0) 
562     phi = phi +2*TMath::Pi();
563   if(phi > 2*TMath::Pi() )
564     phi = phi - 2*TMath::Pi();
565   
566   return phi;
567 }
568 //_____________________________________________________________________
569 Float_t AliFMDAnaParameters::GetEtaFromStrip(UShort_t det, Char_t ring, UShort_t sec, UShort_t strip, Float_t zvtx) const
570 {
571   // AliFMDRing fmdring(ring);
572   // fmdring.Init();
573   Float_t   rad       = GetMaxR(ring)-GetMinR(ring);
574   Float_t   nStrips   = (ring == 'I' ? 512 : 256);
575   Float_t   segment   = rad / nStrips;
576   Float_t   r         = GetMinR(ring) + segment*strip;
577   Float_t   z         = 0;
578   Int_t hybrid = sec / 2;
579   
580   if(det == 1) {
581     if(!(hybrid%2)) z = 320.266; else z = 319.766;
582   }
583   if(det == 2 && ring == 'I' ) {
584     if(!(hybrid%2)) z = 83.666; else z = 83.166;
585   }
586   if(det == 2 && ring == 'O' ) {
587     if(!(hybrid%2)) z = 74.966; else z = 75.466;
588   }
589   if(det == 3 && ring == 'I' ) {
590     if(!(hybrid%2)) z = -63.066; else z = -62.566;
591   }
592   if(det == 3 && ring == 'O' ) {
593     if(!(hybrid%2)) z = -74.966; else z = -75.466;
594   }
595   
596   //std::cout<<det<<"   "<<ring<<"   "<<sec<<"   "<<hybrid<<"    "<<z<<std::endl;
597   
598   // Float_t   r     = TMath::Sqrt(TMath::Power(x,2)+TMath::Power(y,2));
599   Float_t   theta = TMath::ATan2(r,z-zvtx);
600   Float_t   eta   = -1*TMath::Log(TMath::Tan(0.5*theta));
601   
602   return eta;
603 }
604
605 //_____________________________________________________________________
606
607 Bool_t AliFMDAnaParameters::GetVertex(const AliESDEvent* esd, Double_t* vertexXYZ) 
608 {
609   const AliESDVertex* vertex = 0;
610   vertex = esd->GetPrimaryVertexSPD();
611   
612   if(vertex)
613     vertex->GetXYZ(vertexXYZ);
614     
615   //if(vertexXYZ[0] == 0 || vertexXYZ[1] == 0 )
616   //  return kFALSE;
617   
618   if(vertex->GetNContributors() <= 0)
619     return kFALSE;
620   
621   if(vertex->GetZRes() > 0.1 ) 
622     return kFALSE;
623   
624   return vertex->GetStatus();
625   
626 }
627 //____________________________________________________________________
628 Bool_t AliFMDAnaParameters::IsEventTriggered(const AliESDEvent *esd, Trigger trig) {
629
630   Trigger old = fTrigger;
631   fTrigger = trig;
632   Bool_t retval = IsEventTriggered(esd);
633   fTrigger = old;
634   return retval;
635
636 }
637 //____________________________________________________________________
638 Bool_t AliFMDAnaParameters::IsEventTriggered(const AliESDEvent *esd) const {
639   // check if the event was triggered
640   
641   if (fCentralSelection) return kTRUE;
642   ULong64_t triggerMask = esd->GetTriggerMask();
643   
644   TString triggers = esd->GetFiredTriggerClasses();
645   
646   //if(triggers.Contains("CINT1B-ABCE-NOPF-ALL") || triggers.Contains("CINT1B-E-NOPF-ALL")) return kTRUE;
647   //else return kFALSE;
648   //if(triggers.Contains("CINT1B-E-NOPF-ALL"))    return kFALSE;
649   
650   // if(triggers.Contains("CINT1A-ABCE-NOPF-ALL")) return kFALSE;
651   // if(triggers.Contains("CINT1C-ABCE-NOPF-ALL")) return kFALSE;
652   
653   // definitions from p-p.cfg
654   ULong64_t spdFO = (1 << 14);
655   ULong64_t v0left = (1 << 11);
656   ULong64_t v0right = (1 << 12);
657   AliTriggerAnalysis tAna;
658   
659   //REMOVE WHEN FINISHED PLAYING WITH TRIGGERS!
660   //fPhysicsSelection->IsCollisionCandidate(esd);
661   
662   
663   switch (fTrigger) {
664   case kMB1: {
665     if(fRealData) {
666       if( fPhysicsSelection->IsCollisionCandidate(esd))
667         return kTRUE;
668     }
669     else {
670       if (triggerMask & spdFO || ((triggerMask & v0left) || (triggerMask & v0right)))
671         return kTRUE;
672       break;
673     }
674   }
675   case kMB2: { 
676     if (triggerMask & spdFO && ((triggerMask & v0left) || (triggerMask & v0right)))
677       return kTRUE;
678     break;
679   }
680   case kSPDFASTOR: {
681     if (triggerMask & spdFO)
682       return kTRUE;
683     break;
684   }
685   case kNOCTP: {
686     return kTRUE;
687     break;
688   }
689   case kEMPTY: {
690      const AliMultiplicity* testmult = esd->GetMultiplicity();
691     Int_t nTrackLets = testmult->GetNumberOfTracklets();
692     Int_t nClusters  = testmult->GetNumberOfSingleClusters();
693     Int_t fastOR = tAna.SPDFiredChips(esd, 0);
694     
695     const AliESDVertex* vertex = 0;
696     vertex = esd->GetPrimaryVertexSPD();
697     Bool_t vtxstatus = vertex->GetStatus();
698     if(vertex->GetNContributors() <= 0)
699       vtxstatus = kFALSE;
700     
701     if(vertex->GetZRes() > 0.1 ) 
702       vtxstatus = kFALSE;
703     
704     Bool_t v0ABG = tAna.IsOfflineTriggerFired(esd, AliTriggerAnalysis::kV0ABG);
705     
706     Bool_t v0CBG = tAna.IsOfflineTriggerFired(esd, AliTriggerAnalysis::kV0CBG);
707     Bool_t v0A = tAna.IsOfflineTriggerFired(esd, AliTriggerAnalysis::kV0A);
708     Bool_t v0C = tAna.IsOfflineTriggerFired(esd, AliTriggerAnalysis::kV0C);
709     if(triggers.Contains("CINT1A-ABCE-NOPF-ALL") || triggers.Contains("CINT1C-ABCE-NOPF-ALL")) 
710       return kTRUE;
711     break;
712   }
713     
714   }//switch
715   
716   return kFALSE;
717 }
718     
719 //____________________________________________________________________
720 Float_t 
721 AliFMDAnaParameters::GetStripLength(Char_t ring, UShort_t strip)  
722 {
723   //AliFMDRing fmdring(ring);
724   // fmdring.Init();
725   
726   Float_t rad        = GetMaxR(ring)-GetMinR(ring);
727   Float_t   nStrips   = (ring == 'I' ? 512 : 256);
728   Float_t segment    = rad / nStrips;
729   
730   //TVector2* corner1  = fmdring.GetVertex(2);  
731   // TVector2* corner2  = fmdring.GetVertex(3);
732   
733   SetCorners(ring);
734   /*
735   std::cout<<GetMaxR(ring)<<"   "<<fmdring.GetMaxR()<<std::endl;
736   std::cout<<GetMinR(ring)<<"   "<<fmdring.GetMinR()<<std::endl;
737   std::cout<<corner1->X()<<"   "<<fCorner1.X()<<std::endl;
738   std::cout<<corner2->X()<<"   "<<fCorner2.X()<<std::endl;
739   std::cout<<corner1->Y()<<"   "<<fCorner1.Y()<<std::endl;
740   std::cout<<corner2->Y()<<"   "<<fCorner2.Y()<<std::endl;*/
741   Float_t slope      = (fCorner1.Y() - fCorner2.Y()) / (fCorner1.X() - fCorner2.X());
742   Float_t constant   = (fCorner2.Y()*fCorner1.X()-(fCorner2.X()*fCorner1.Y())) / (fCorner1.X() - fCorner2.X());
743   Float_t radius     = GetMinR(ring) + strip*segment;
744   
745   Float_t d          = TMath::Power(TMath::Abs(radius*slope),2) + TMath::Power(radius,2) - TMath::Power(constant,2);
746   
747   Float_t arclength  = GetBaseStripLength(ring,strip);
748   if(d>0) {
749     
750     Float_t x        = (-1*TMath::Sqrt(d) -slope*constant) / (1+TMath::Power(slope,2));
751     Float_t y        = slope*x + constant;
752     Float_t theta    = TMath::ATan2(x,y);
753     
754     if(x < fCorner1.X() && y > fCorner1.Y()) {
755       arclength = radius*theta;                        //One sector since theta is by definition half-hybrid
756       
757     }
758     
759   }
760   
761   return arclength;
762   
763   
764 }
765 //____________________________________________________________________
766 Float_t 
767 AliFMDAnaParameters::GetBaseStripLength(Char_t ring, UShort_t strip)  
768 {  
769   // AliFMDRing fmdring(ring);
770   // fmdring.Init();
771   Float_t rad             = GetMaxR(ring)-GetMinR(ring);
772   Float_t nStrips         = (ring == 'I' ? 512 : 256);
773   Float_t nSec            = (ring == 'I' ? 20 : 40);
774   Float_t segment         = rad / nStrips;
775   Float_t basearc         = 2*TMath::Pi() / (0.5*nSec); // One hybrid: 36 degrees inner, 18 outer
776   Float_t radius          = GetMinR(ring) + strip*segment;
777   Float_t basearclength   = 0.5*basearc * radius;                // One sector   
778   
779   return basearclength;
780 }
781 //____________________________________________________________________
782 //
783 // EOF
784 //