]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/HESDfromKin.C
minors.
[u/mrichter/AliRoot.git] / HMPID / HESDfromKin.C
1 AliRunLoader *gAL=0; 
2 Int_t gEvt=0; Int_t gMaxEvt=0;
3 TObjArray *pNmean,*pQthre;
4 TTree *gEsdTr;
5 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
6 void HESDfromKin(const char *name="default")
7 {//simulate ESD from kinematics
8
9   if(gSystem->IsFileInIncludePath("galice.root")){// tries to open session
10     if(gAlice) delete gAlice;                                               //in case we execute this in aliroot delete default AliRun object 
11     gAL=AliRunLoader::Open();                                                                    //try to open galice.root from current dir 
12     gAL->LoadgAlice();                                                                           //take new AliRun object from galice.root   
13     if(gAL->LoadHeader()) return;
14     if(gAL->LoadKinematics()) return;
15
16     AliLoader *pHL=gAL->GetDetectorLoader("HMPID");
17     pHL->LoadRecPoints();
18     AliESDEvent *pEsd = new AliESDEvent();   
19     TFile *pEsdFl=TFile::Open("AliESDs.root","recreate"); 
20     gEsdTr=new TTree("esdTree","Sim ESD from kinematics"); 
21     pEsd->CreateStdContent();    pEsd->WriteToTree(gEsdTr);  //clm: new ESD write schema: see Task Force meeting 20th June, 2007
22     gEsdTr->GetUserInfo()->Add(pEsd);                        //clm: TList has to be created for ReadFromTree method -- this was not needed by the old ESD
23  
24        
25   }  else return;  
26
27   if(!OpenCalib()) {Printf("Problems in OpenCalib!Bye.");return;}
28     
29   TString ttl=name;
30   Bool_t htaCheck=ttl.Contains("HTA");
31 //  if(!htaCheck) SimEsd(pHL,pEsd); else SimEsdHidden(pHL,pEsd);
32   SimEsd(pHL,pEsd,htaCheck);
33   
34   pEsdFl->cd();
35   pEsdFl->Write();pEsdFl->Close();        
36 }
37 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
38 void SimEsd(AliLoader *pHL,AliESDEvent *pEsd,Bool_t htaCheck)
39 {
40   if(htaCheck) {
41     TFile *fout = new TFile("HTA.root","recreate");
42     TH1F *hdC   = new TH1F("dC"  ,";delta Cerenkov (rad)",100,-0.2,0.2);
43     TH1F *hCer  = new TH1F("Cer" ,"Theta Cerenkov (rad)",250,0.,0.75);
44     TH2F *htvsp = new TH2F("tvsp",";momentum (GeV/c);theta Cerenkov (rad)",100,0.,5.,1000,0.,0.75);
45     TH1F *hdth  = new TH1F("dth" ,";Delta theta Trk (mrad)",100,-250,250);
46     TH1F *hdph  = new TH1F("dph" ,";Delta phi Trk (mrad)",100,-500,500);
47     Double_t rd=TMath::RadToDeg();
48     Printf("----------------------------------------------");
49     Printf("| SimHTA:Utility to embed ESD from kinematics|");
50     Printf("|     with  Hidden Track Algorithm (HTA)     |");
51     Printf("----------------------------------------------");
52   } else {
53     Printf("-----------------------------------------------");
54     Printf("| SimESD: Utility to embed ESD from kinematics|");
55     Printf("-----------------------------------------------");
56 }
57   AliHMPIDTracker::SetFieldMap(gAL->GetAliRun()->Field(),kTRUE);
58   AliHMPID *pH=(AliHMPID*)gAL->GetAliRun()->GetDetector("HMPID");
59   Int_t iNevt=gAL->GetNumberOfEvents();
60   pEsd->SetMagneticField(AliHMPIDTracker::GetBz());
61   for(Int_t iEvt=0;iEvt<iNevt;iEvt++){//events loop
62     gAL->GetEvent(iEvt);    
63     pHL->TreeR()->GetEntry(0);
64     AliStack *pStack=gAL->Stack();
65     Int_t nTrkHMPID=0;
66     for(Int_t i=0;i<pStack->GetNtrack();i++){
67       if(!pStack->IsPhysicalPrimary(i)) continue;
68       TParticle *pTrack=pStack->Particle(i); 
69       if(pTrack->GetPDG()->Charge()==0) continue;
70       AliESDtrack trk(pTrack); 
71       Float_t xPc,yPc,xRa,yRa,thRa,phRa;
72       Int_t iCh=AliHMPIDTracker::IntTrkCha(&trk,xPc,yPc,xRa,yRa,thRa,phRa);         //get chamber intersected by this track 
73       if(iCh<0) {
74         trk.SetHMPIDtrk(0,0,0,0);                                                                //no intersection found
75         trk.SetHMPIDcluIdx   (99,99999);                                                         //chamber not found, mip not yet considered
76         trk.SetHMPIDsignal(AliHMPIDRecon::kNotPerformed);                                        //ring reconstruction not yet performed
77         continue;                                                           //no intersection at all, go after next track
78       }
79       nTrkHMPID++;
80       trk.SetHMPIDcluIdx   (iCh,99999);                                                          //chamber not found, mip not yet considered
81       trk.SetHMPIDtrk(xRa,yRa,thRa,phRa);                                                        //store initial infos
82       pEsd->AddTrack(&trk);
83     }// track loop
84     
85     if(!(iEvt%50)) Printf("Number of events processed: %i with tracks %i in HMPID",iEvt,nTrkHMPID);
86     if(!htaCheck) AliHMPIDTracker::Recon(pEsd,pH->CluLst(),pNmean,pQthre);
87 // Perform PID
88     
89     
90     for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){                                       //ESD tracks loop
91       
92       AliESDtrack *pTrk = pEsd->GetTrack(iTrk);                                                    //get reconstructed track    
93       
94       Float_t radX,radY,thetaTrk,phiTrk;
95       pTrk->GetHMPIDtrk(radX,radY,thetaTrk,phiTrk);
96       
97       if(htaCheck) {
98           Int_t iCh = pTrk->GetHMPIDcluIdx();                                                                  //chamber not found, mip not yet considered
99           iCh/=1000000;
100         Printf("simulated track theta %f phi %f",thetaTrk*rd,phiTrk*rd);
101         TObjArray *pClus = pH->CluLst();
102
103         if(AliHMPIDTracker::ReconHiddenTrk(pEsd,pClus,pNmean,pQthre)!=0) continue;
104
105         Double_t thetaCerSim = TMath::ACos(pTrack->Energy()/(1.292*pTrack->P()));
106         Printf("-----------------------------------------------------------");
107         Printf(" theta Cerenkov simulated     %f",thetaCerSim);
108         Printf(" theta Cerenkov reconstructed %f",pTrk->GetHMPIDsignal());
109         Float_t thetaFit,phiFit;
110         pTrk->GetHMPIDtrk(radX,radY,thetaFit,phiFit);
111 //        Printf("reconstr. track theta %f phi %f",thetaFit*rd,phiFit*rd);
112
113         // fill useful histos
114         hdC->Fill(pTrk->GetHMPIDsignal()-thetaCerSim);
115         hCer->Fill(pTrk->GetHMPIDsignal());
116         htvsp->Fill(pTrk->P(),pTrk->GetHMPIDsignal());
117         hdth->Fill((thetaFit-thetaTrk)*1000);
118         hdph->Fill((phiFit-phiTrk)*1000);
119       }
120       
121       AliHMPIDPid pID;
122       Double_t prob[5];
123       pID.FindPid(pTrk,5,prob);
124       pTrk->SetHMPIDpid(prob);
125 //      Printf(" Prob e- %6.2f mu %6.2f pi %6.2f k %6.2f p %6.2f",prob[0]*100,prob[1]*100,prob[2]*100,prob[3]*100,prob[4]*100);
126     }
127     
128     gEsdTr->Fill();
129     pEsd->Reset();
130     
131   }// event loop
132   Printf("Events processed %i",iEvt);
133   if(htaCheck) {
134     fout->Write();
135     fout->Close();
136     delete fout;
137   }
138   gAL->UnloadHeader();  gAL->UnloadKinematics();
139 }//Esd()
140 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
141 void SimEsdHidden(AliLoader *pHL,AliESDEvent *pEsd)
142 {
143   TFile *fout = new TFile("HTA.root","recreate");
144   TH1F *hdC   = new TH1F("dC"  ,";delta Cerenkov (rad)",100,-0.2,0.2);
145   TH1F *hCer  = new TH1F("Cer" ,"Theta Cerenkov (rad)",250,0.,0.75);
146   TH2F *htvsp = new TH2F("tvsp",";momentum (GeV/c);theta Cerenkov (rad)",100,0.,5.,1000,0.,0.75);
147   TH1F *hdth  = new TH1F("dth" ,";Delta theta Trk (mrad)",100,-250,250);
148   TH1F *hdph  = new TH1F("dph" ,";Delta phi Trk (mrad)",100,-500,500);
149   Double_t rd=TMath::RadToDeg();
150   Printf("----------------------------------------------");
151   Printf("| SimHTA:Utility to embed ESD from kinematics|");
152   Printf("|     with  Hidden Track Algorithm (HTA)     |");
153   Printf("----------------------------------------------");
154   AliHMPIDTracker::SetFieldMap(gAL->GetAliRun()->Field(),kTRUE);
155   AliHMPID *pH=(AliHMPID*)gAL->GetAliRun()->GetDetector("HMPID");
156   Int_t iNevt=gAL->GetNumberOfEvents();
157   pEsd->SetMagneticField(AliHMPIDTracker::GetBz());
158   
159   Printf("Number of events to process: %i",iNevt);
160   
161   for(Int_t iEvt=0;iEvt<iNevt;iEvt++){//events loop
162     if(!(iEvt%1)) Printf("Events processed %i",iEvt);
163     gAL->GetEvent(iEvt);    
164     pHL->TreeR()->GetEntry(0);
165     AliStack *pStack=gAL->Stack();
166     
167     for(Int_t i=0;i<pStack->GetNtrack();i++){
168       
169       (!pStack->IsPhysicalPrimary(i)) continue;
170       TParticle *pTrack=pStack->Particle(i); 
171       
172       if(pTrack->GetPDG()->Charge()==0) continue;
173       //find the chamber that intersects HMPID
174       pTrack->Print();
175       AliESDtrack trk(pTrack);
176       Float_t xPc,yPc,xRa,yRa,thRa,phRa;
177       Int_t iCh=AliHMPIDTracker::IntTrkCha(&trk,xPc,yPc,xRa,yRa,thRa,phRa);         //get chamber intersected by this track
178       if(iCh<0) continue;                                                           //no intersection at all, go after next track
179       trk.SetHMPIDtrk(xRa,yRa,thRa,phRa);                                                        //store initial infos
180       Float_t radX,radY,thetaTrk,phiTrk;
181       trk.GetHMPIDtrk(radX,radY,thetaTrk,phiTrk);
182 //      Printf("simulated track theta %f phi %f",thetaTrk*rd,phiTrk*rd);
183       TObjArray *pClus = pH->CluLst();
184       
185       if(AliHMPIDTracker::ReconHiddenTrk(iCh,0,&trk,(TClonesArray *)pClus->At(iCh),pNmean,pQthre)!=0) continue;
186       
187       pEsd->AddTrack(&trk);
188       Double_t thetaCerSim = TMath::ACos(pTrack->Energy()/(1.292*pTrack->P()));
189 //      Printf(" theta Cerenkov simulated     %f",thetaCerSim);
190 //      Printf(" theta Cerenkov reconstructed %f",trk.GetHMPIDsignal());
191       Float_t thetaFit,phiFit;
192       trk.GetHMPIDtrk(radX,radY,thetaFit,phiFit);
193 //      Printf("reconstr. track theta %f phi %f",thetaFit*rd,phiFit*rd);
194       
195       // fill useful histos
196       hdC->Fill(trk.GetHMPIDsignal()-thetaCerSim);
197       hCer->Fill(trk.GetHMPIDsignal());
198       htvsp->Fill(pTrack->P(),trk.GetHMPIDsignal());
199       hdth->Fill((thetaFit-thetaTrk)*1000);
200       hdph->Fill((phiFit-phiTrk)*1000);
201     }// track loop
202     pEsd->SetMagneticField(AliHMPIDTracker::GetBz());
203     gEsdTr->Fill();
204     pEsd->Reset();
205   }// event loop
206   Printf("Events processed %i",iEvt);
207   fout->Write();
208   fout->Close();
209   delete fout;
210   gAL->UnloadHeader();  gAL->UnloadKinematics();
211 }//EsdHidden()
212 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
213 Bool_t OpenCalib()
214 {
215   AliCDBManager* pCDB = AliCDBManager::Instance();
216   pCDB->SetDefaultStorage("local://$HOME");
217   AliCDBEntry *pQthreEnt=pCDB->Get("HMPID/Calib/Qthre",0);
218   AliCDBEntry *pNmeanEnt=pCDB->Get("HMPID/Calib/Nmean",0);
219   
220   if(!pQthreEnt || !pNmeanEnt) return kFALSE;
221   
222   pNmean=(TObjArray*)pNmeanEnt->GetObject(); 
223   pQthre=(TObjArray*)pQthreEnt->GetObject(); 
224
225   if(!pQthre || !pNmean) return kFALSE;  
226   return kTRUE;
227 }
228 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++