]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/HESDfromKin.C
Comments corrected.
[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   
33   pEsdFl->cd();
34   pEsdFl->Write();pEsdFl->Close();        
35 }
36 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
37 void SimEsd(AliLoader *pHL,AliESDEvent *pEsd)
38 {
39   Printf("-----------------------------------------------");
40   Printf("| SimESD: Utility to embed ESD from kinematics|");
41   Printf("-----------------------------------------------");
42   AliHMPIDTracker::SetFieldMap(gAL->GetAliRun()->Field(),kTRUE);
43   AliHMPID *pH=(AliHMPID*)gAL->GetAliRun()->GetDetector("HMPID");
44   Int_t mtid=-1;
45   Int_t iNevt=gAL->GetNumberOfEvents();
46   pEsd->SetMagneticField(AliHMPIDTracker::GetBz());
47   Printf("Number of events to process: %i",iNevt);
48   for(Int_t iEvt=0;iEvt<iNevt;iEvt++){//events loop
49     if(!(iEvt%50)) Printf("Events processed %i",iEvt);
50     gAL->GetEvent(iEvt);    
51     pHL->TreeR()->GetEntry(0);
52     AliStack *pStack=gAL->Stack();
53     for(Int_t i=0;i<pStack->GetNtrack();i++){
54       TParticle *pTrack=pStack->Particle(i); 
55       mtid=pTrack->GetFirstMother();
56       if(mtid>=0) continue; // only primaries
57       if(pTrack->GetPDG()->Charge()==0) continue;
58       AliESDtrack trk(pTrack); 
59       Float_t xPc,yPc;
60       Int_t iCh=AliHMPIDTracker::IntTrkCha(&trk,xPc,yPc);                           //get chamber intersected by this track 
61       if(iCh<0) continue;                                                           //no intersection at all, go after next track
62 //      Printf(" particle analyzed %d with mtid %d is %s hitting chamber %d",i,mtid,pTrack->GetName(),iCh);
63       pEsd->AddTrack(&trk);
64       AliHMPIDTracker::Recon(pEsd,pH->CluLst(),pNmean,pQthre);
65     }// track loop
66     gEsdTr->Fill();
67     pEsd->Reset();
68   }// event loop
69   Printf("Events processed %i",iEvt);
70   gAL->UnloadHeader();  gAL->UnloadKinematics();
71 }//Esd()
72 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
73 void SimEsdHidden(AliLoader *pHL,AliESDEvent *pEsd)
74 {
75   TFile *fout = new TFile("HTA.root","recreate");
76   TH1F *hdC   = new TH1F("dC"  ,";delta Cerenkov (rad)",100,-0.2,0.2);
77   TH1F *hCer  = new TH1F("Cer" ,"Theta Cerenkov (rad)",250,0.,0.75);
78   TH2F *htvsp = new TH2F("tvsp",";momentum (GeV/c);theta Cerenkov (rad)",100,0.,5.,1000,0.,0.75);
79   TH1F *hdth  = new TH1F("dth" ,";Delta theta Trk (mrad)",100,-250,250);
80   TH1F *hdph  = new TH1F("dph" ,";Delta phi Trk (mrad)",100,-500,500);
81   Double_t rd=TMath::RadToDeg();
82   Printf("----------------------------------------------");
83   Printf("| SimHTA:Utility to embed ESD from kinematics|");
84   Printf("|     with  Hidden Track Algorithm (HTA)     |");
85   Printf("----------------------------------------------");
86   AliHMPIDTracker::SetFieldMap(gAL->GetAliRun()->Field(),kTRUE);
87   AliHMPID *pH=(AliHMPID*)gAL->GetAliRun()->GetDetector("HMPID");
88   Int_t mtid=-1;
89   Int_t iNevt=gAL->GetNumberOfEvents();
90   
91   Printf("Number of events to process: %i",iNevt);
92   
93   for(Int_t iEvt=0;iEvt<iNevt;iEvt++){//events loop
94     if(!(iEvt%50)) Printf("Events processed %i",iEvt);
95     gAL->GetEvent(iEvt);    
96     pHL->TreeR()->GetEntry(0);
97     AliStack *pStack=gAL->Stack();
98     
99     for(Int_t i=0;i<pStack->GetNtrack();i++){
100       
101       TParticle *pTrack=pStack->Particle(i); 
102       mtid=pTrack->GetFirstMother();
103       
104       if(mtid>=0) continue; // only primaries
105       
106       //find the chamber that intersects HMPID
107       AliESDtrack trk(pTrack);
108       Float_t xPc,yPc;
109       Int_t iCh=AliHMPIDTracker::IntTrkCha(&trk,xPc,yPc);                           //get chamber intersected by this track 
110       if(iCh<0) continue;                                                           //no intersection at all, go after next track
111       Float_t radX,radY,thetaTrk,phiTrk;
112       trk.GetHMPIDtrk(radX,radY,thetaTrk,phiTrk);
113 //      Printf("simulated track theta %f phi %f",thetaTrk*rd,phiTrk*rd);
114       TObjArray *pClus = pH->CluLst();
115       
116       if(AliHMPIDTracker::ReconHiddenTrk(iCh,&trk,(TClonesArray *)pClus->At(iCh),pNmean,pQthre)!=0) continue;
117       
118       pEsd->AddTrack(&trk);
119       Double_t thetaCerSim = TMath::ACos(pTrack->Energy()/(1.292*pTrack->P()));
120 //      Printf(" theta Cerenkov simulated     %f",thetaCerSim);
121 //      Printf(" theta Cerenkov reconstructed %f",trk.GetHMPIDsignal());
122       Float_t thetaFit,phiFit;
123       trk.GetHMPIDtrk(radX,radY,thetaFit,phiFit);
124 //      Printf("reconstr. track theta %f phi %f",thetaFit*rd,phiFit*rd);
125       
126       // fill useful histos
127       hdC->Fill(trk.GetHMPIDsignal()-thetaCerSim);
128       hCer->Fill(trk.GetHMPIDsignal());
129       htvsp->Fill(pTrack->P(),trk.GetHMPIDsignal());
130       hdth->Fill((thetaFit-thetaTrk)*1000);
131       hdph->Fill((phiFit-phiTrk)*1000);
132     }// track loop
133     pEsd->SetMagneticField(AliHMPIDTracker::GetBz());
134     gEsdTr->Fill();
135     pEsd->Reset();
136   }// event loop
137   Printf("Events processed %i",iEvt);
138   fout->Write();
139   fout->Close();
140   delete fout;
141   gAL->UnloadHeader();  gAL->UnloadKinematics();
142 }//EsdHidden()
143 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
144 Bool_t OpenCalib()
145 {
146   AliCDBManager* pCDB = AliCDBManager::Instance();
147   pCDB->SetDefaultStorage("local://$HOME");
148   AliCDBEntry *pQthreEnt=pCDB->Get("HMPID/Calib/Qthre",0);
149   AliCDBEntry *pNmeanEnt=pCDB->Get("HMPID/Calib/Nmean",0);
150   
151   if(!pQthreEnt || !pNmeanEnt) return kFALSE;
152   
153   pNmean=(TObjArray*)pNmeanEnt->GetObject(); 
154   pQthre=(TObjArray*)pQthreEnt->GetObject(); 
155
156   if(!pQthre || !pNmean) return kFALSE;  
157   return kTRUE;
158 }
159 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++