]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/HESDfromKin.C
Comment...uncomment...
[u/mrichter/AliRoot.git] / HMPID / HESDfromKin.C
CommitLineData
d8770a0f 1AliRunLoader *gAL=0;
2Int_t gEvt=0; Int_t gMaxEvt=0;
3df18479 3TObjArray *pNmean,*pQthre;
d8770a0f 4TTree *gEsdTr;
5//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
e4290ede 6void HESDfromKin(const char *name="default")
d8770a0f 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();
6e06db1d 18 AliESDEvent *pEsd = new AliESDEvent();
d8770a0f 19 TFile *pEsdFl=TFile::Open("AliESDs.root","recreate");
20 gEsdTr=new TTree("esdTree","Sim ESD from kinematics");
e4290ede 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
d8770a0f 25 } else return;
26
3df18479 27 if(!OpenCalib()) {Printf("Problems in OpenCalib!Bye.");return;}
d8770a0f 28
c0d7adf8 29 TString ttl=name;
30 Bool_t htaCheck=ttl.Contains("HTA");
9785d5fb 31// if(!htaCheck) SimEsd(pHL,pEsd); else SimEsdHidden(pHL,pEsd);
32 SimEsd(pHL,pEsd,htaCheck);
e4290ede 33
34 pEsdFl->cd();
d8770a0f 35 pEsdFl->Write();pEsdFl->Close();
36}
37//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
9785d5fb 38void SimEsd(AliLoader *pHL,AliESDEvent *pEsd,Bool_t htaCheck)
d8770a0f 39{
9785d5fb 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}
d8770a0f 57 AliHMPIDTracker::SetFieldMap(gAL->GetAliRun()->Field(),kTRUE);
e56b695f 58 AliHMPIDTracker pTracker;
d8770a0f 59 AliHMPID *pH=(AliHMPID*)gAL->GetAliRun()->GetDetector("HMPID");
d8770a0f 60 Int_t iNevt=gAL->GetNumberOfEvents();
1b46ef29 61 pEsd->SetMagneticField(AliHMPIDTracker::GetBz());
d8770a0f 62 for(Int_t iEvt=0;iEvt<iNevt;iEvt++){//events loop
d8770a0f 63 gAL->GetEvent(iEvt);
64 pHL->TreeR()->GetEntry(0);
65 AliStack *pStack=gAL->Stack();
9785d5fb 66 Int_t nTrkHMPID=0;
2b919fd5 67
d8770a0f 68 for(Int_t i=0;i<pStack->GetNtrack();i++){
9785d5fb 69 if(!pStack->IsPhysicalPrimary(i)) continue;
d8770a0f 70 TParticle *pTrack=pStack->Particle(i);
1b46ef29 71 if(pTrack->GetPDG()->Charge()==0) continue;
e4290ede 72 AliESDtrack trk(pTrack);
39cd22e6 73 Float_t xPc,yPc,xRa,yRa,thRa,phRa;
e56b695f 74 Int_t iCh=pTracker.IntTrkCha(&trk,xPc,yPc,xRa,yRa,thRa,phRa); //get chamber intersected by this track
39cd22e6 75 if(iCh<0) {
76 trk.SetHMPIDtrk(0,0,0,0); //no intersection found
77 trk.SetHMPIDcluIdx (99,99999); //chamber not found, mip not yet considered
78 trk.SetHMPIDsignal(AliHMPIDRecon::kNotPerformed); //ring reconstruction not yet performed
79 continue; //no intersection at all, go after next track
80 }
9785d5fb 81 nTrkHMPID++;
39cd22e6 82 trk.SetHMPIDcluIdx (iCh,99999); //chamber not found, mip not yet considered
e56b695f 83
84 if(phRa<0) phRa += TMath::TwoPi(); // to be verified
85
39cd22e6 86 trk.SetHMPIDtrk(xRa,yRa,thRa,phRa); //store initial infos
d8770a0f 87 pEsd->AddTrack(&trk);
9785d5fb 88
6574baa2 89 Int_t status;
90 if(!htaCheck) status = pTracker.Recon (pEsd,pH->CluLst(),pNmean,pQthre);
91 else status = pTracker.ReconHiddenTrk(pEsd,pH->CluLst(),pNmean,pQthre);
92
2b919fd5 93// Printf("status %i",status);
6574baa2 94 if(status !=0) continue;
9785d5fb 95
ac66a50d 96
6574baa2 97 }// track loop
2b919fd5 98
99 if(!(iEvt%50)) Printf("Number of events processed: %i with tracks %i in HMPID",iEvt,nTrkHMPID);
100// Printf("Number of events processed: %i with tracks %i in HMPID",iEvt,nTrkHMPID);
101
d8770a0f 102 gEsdTr->Fill();
103 pEsd->Reset();
104 }// event loop
105 Printf("Events processed %i",iEvt);
204299b8 106 if(htaCheck) {
107 fout->Write();
108 fout->Close();
109 delete fout;
110 }
d8770a0f 111 gAL->UnloadHeader(); gAL->UnloadKinematics();
e56b695f 112}//SimEsd()
d8770a0f 113//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3df18479 114Bool_t OpenCalib()
d8770a0f 115{
116 AliCDBManager* pCDB = AliCDBManager::Instance();
117 pCDB->SetDefaultStorage("local://$HOME");
118 AliCDBEntry *pQthreEnt=pCDB->Get("HMPID/Calib/Qthre",0);
119 AliCDBEntry *pNmeanEnt=pCDB->Get("HMPID/Calib/Nmean",0);
120
3df18479 121 if(!pQthreEnt || !pNmeanEnt) return kFALSE;
d8770a0f 122
123 pNmean=(TObjArray*)pNmeanEnt->GetObject();
3df18479 124 pQthre=(TObjArray*)pQthreEnt->GetObject();
125
126 if(!pQthre || !pNmean) return kFALSE;
127 return kTRUE;
d8770a0f 128}
129//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++