ROOT logo
AliRunLoader *gAL=0; 
Int_t gEvt=0; Int_t gMaxEvt=0;
TObjArray *pNmean,*pQthre;
TTree *gEsdTr;
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void HESDfromKin(const char *name="default")
{//simulate ESD from kinematics

  if(gSystem->IsFileInIncludePath("galice.root")){// tries to open session
    if(gAlice) delete gAlice;                                               //in case we execute this in aliroot delete default AliRun object 
    gAL=AliRunLoader::Open();                                                                    //try to open galice.root from current dir 
    gAL->LoadgAlice();                                                                           //take new AliRun object from galice.root   
    if(gAL->LoadHeader()) return;
    if(gAL->LoadKinematics()) return;

    AliLoader *pHL=gAL->GetDetectorLoader("HMPID");
    pHL->LoadRecPoints();
    AliESDEvent *pEsd = new AliESDEvent();   
    TFile *pEsdFl=TFile::Open("AliESDs.root","recreate"); 
    gEsdTr=new TTree("esdTree","Sim ESD from kinematics"); 
    pEsd->CreateStdContent();    pEsd->WriteToTree(gEsdTr);  //clm: new ESD write schema: see Task Force meeting 20th June, 2007
    gEsdTr->GetUserInfo()->Add(pEsd);                        //clm: TList has to be created for ReadFromTree method -- this was not needed by the old ESD
 
       
  }  else return;  

  if(!OpenCalib()) {Printf("Problems in OpenCalib!Bye.");return;}
    
  TString ttl=name;
  Bool_t htaCheck=ttl.Contains("HTA");
//  if(!htaCheck) SimEsd(pHL,pEsd); else SimEsdHidden(pHL,pEsd);
  SimEsd(pHL,pEsd,htaCheck);
  
  pEsdFl->cd();
  pEsdFl->Write();pEsdFl->Close();        
}
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
void SimEsd(AliLoader *pHL,AliESDEvent *pEsd,Bool_t htaCheck)
{
  if(htaCheck) {
    TFile *fout = new TFile("HTA.root","recreate");
    TH1F *hdC   = new TH1F("dC"  ,";delta Cerenkov (rad)",100,-0.2,0.2);
    TH1F *hCer  = new TH1F("Cer" ,"Theta Cerenkov (rad)",250,0.,0.75);
    TH2F *htvsp = new TH2F("tvsp",";momentum (GeV/c);theta Cerenkov (rad)",100,0.,5.,1000,0.,0.75);
    TH1F *hdth  = new TH1F("dth" ,";Delta theta Trk (mrad)",100,-250,250);
    TH1F *hdph  = new TH1F("dph" ,";Delta phi Trk (mrad)",100,-500,500);
    Double_t rd=TMath::RadToDeg();
    Printf("----------------------------------------------");
    Printf("| SimHTA:Utility to embed ESD from kinematics|");
    Printf("|     with  Hidden Track Algorithm (HTA)     |");
    Printf("----------------------------------------------");
  } else {
    Printf("-----------------------------------------------");
    Printf("| SimESD: Utility to embed ESD from kinematics|");
    Printf("-----------------------------------------------");
}

  AliGRPManager *grpMan = new AliGRPManager();
  
  grpMan->ReadGRPEntry();
  grpMan->SetMagField();

  AliHMPIDTracker pTracker;
  AliHMPID *pH=(AliHMPID*)gAL->GetAliRun()->GetDetector("HMPID");
  Int_t iNevt=gAL->GetNumberOfEvents();
  pEsd->SetMagneticField(AliHMPIDTracker::GetBz());
  for(Int_t iEvt=0;iEvt<iNevt;iEvt++){//events loop
    gAL->GetEvent(iEvt);    
    pHL->TreeR()->GetEntry(0);
    AliStack *pStack=gAL->Stack();
    Int_t nTrkHMPID=0;

    for(Int_t i=0;i<pStack->GetNtrack();i++){
      if(!pStack->IsPhysicalPrimary(i)) continue;
      TParticle *pTrack=pStack->Particle(i); 
      if(pTrack->GetPDG()->Charge()==0) continue;
//      Printf("track n. %i",i);
      AliESDtrack trk(pTrack); 
      Float_t xPc,yPc,xRa,yRa,thRa,phRa;
      Int_t iCh=pTracker.IntTrkCha(&trk,xPc,yPc,xRa,yRa,thRa,phRa);         //get chamber intersected by this track 
      if(iCh<0) {
        trk.SetHMPIDtrk(0,0,0,0);                                                                //no intersection found
        trk.SetHMPIDcluIdx   (99,99999);                                                         //chamber not found, mip not yet considered
        trk.SetHMPIDsignal(AliHMPIDRecon::kNotPerformed);                                        //ring reconstruction not yet performed
        continue;                                                           //no intersection at all, go after next track
      }
      nTrkHMPID++;
      trk.SetHMPIDcluIdx   (iCh,99999);                                                          //chamber not found, mip not yet considered
      
      if(phRa<0) phRa += TMath::TwoPi(); // to be verified
      
      trk.SetHMPIDtrk(xPc,yPc,thRa,phRa);                                                        //store initial infos
      trk.SetLabel(i);
      pEsd->AddTrack(&trk);
    
      Int_t status;
      if(!htaCheck) status = pTracker.ReconFromKin  (pEsd,pH->CluLst(),pNmean,pQthre);
      else          status = pTracker.ReconHiddenTrk(pEsd,pH->CluLst(),pNmean,pQthre);

//      Printf("status %i",status);
      if(status !=0) continue;

    
    }// track loop
    
    if(!(iEvt%50)) Printf("Number of events processed: %i with tracks %i in HMPID",iEvt,nTrkHMPID);
//      Printf("Number of events processed: %i with tracks %i in HMPID",iEvt,nTrkHMPID);
    
    gEsdTr->Fill();
    pEsd->Reset();
  }// event loop
  Printf("Events processed %i",iEvt);
  if(htaCheck) {
    fout->Write();
    fout->Close();
    delete fout;
  }
  gAL->UnloadHeader();  gAL->UnloadKinematics();
}//SimEsd()
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Bool_t OpenCalib()
{
  AliCDBManager* pCDB = AliCDBManager::Instance();
  pCDB->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
  pCDB->SetSpecificStorage("GRP/GRP/Data", Form("local://%s",gSystem->pwd()));
  pCDB->SetRun(0);
  AliCDBEntry *pQthreEnt=pCDB->Get("HMPID/Calib/Qthre",0);
  AliCDBEntry *pNmeanEnt=pCDB->Get("HMPID/Calib/Nmean",0);
  
  if(!pQthreEnt || !pNmeanEnt) return kFALSE;
  
  pNmean=(TObjArray*)pNmeanEnt->GetObject(); 
  pQthre=(TObjArray*)pQthreEnt->GetObject(); 

  if(!pQthre || !pNmean) return kFALSE;  
  return kTRUE;
}
 HESDfromKin.C:1
 HESDfromKin.C:2
 HESDfromKin.C:3
 HESDfromKin.C:4
 HESDfromKin.C:5
 HESDfromKin.C:6
 HESDfromKin.C:7
 HESDfromKin.C:8
 HESDfromKin.C:9
 HESDfromKin.C:10
 HESDfromKin.C:11
 HESDfromKin.C:12
 HESDfromKin.C:13
 HESDfromKin.C:14
 HESDfromKin.C:15
 HESDfromKin.C:16
 HESDfromKin.C:17
 HESDfromKin.C:18
 HESDfromKin.C:19
 HESDfromKin.C:20
 HESDfromKin.C:21
 HESDfromKin.C:22
 HESDfromKin.C:23
 HESDfromKin.C:24
 HESDfromKin.C:25
 HESDfromKin.C:26
 HESDfromKin.C:27
 HESDfromKin.C:28
 HESDfromKin.C:29
 HESDfromKin.C:30
 HESDfromKin.C:31
 HESDfromKin.C:32
 HESDfromKin.C:33
 HESDfromKin.C:34
 HESDfromKin.C:35
 HESDfromKin.C:36
 HESDfromKin.C:37
 HESDfromKin.C:38
 HESDfromKin.C:39
 HESDfromKin.C:40
 HESDfromKin.C:41
 HESDfromKin.C:42
 HESDfromKin.C:43
 HESDfromKin.C:44
 HESDfromKin.C:45
 HESDfromKin.C:46
 HESDfromKin.C:47
 HESDfromKin.C:48
 HESDfromKin.C:49
 HESDfromKin.C:50
 HESDfromKin.C:51
 HESDfromKin.C:52
 HESDfromKin.C:53
 HESDfromKin.C:54
 HESDfromKin.C:55
 HESDfromKin.C:56
 HESDfromKin.C:57
 HESDfromKin.C:58
 HESDfromKin.C:59
 HESDfromKin.C:60
 HESDfromKin.C:61
 HESDfromKin.C:62
 HESDfromKin.C:63
 HESDfromKin.C:64
 HESDfromKin.C:65
 HESDfromKin.C:66
 HESDfromKin.C:67
 HESDfromKin.C:68
 HESDfromKin.C:69
 HESDfromKin.C:70
 HESDfromKin.C:71
 HESDfromKin.C:72
 HESDfromKin.C:73
 HESDfromKin.C:74
 HESDfromKin.C:75
 HESDfromKin.C:76
 HESDfromKin.C:77
 HESDfromKin.C:78
 HESDfromKin.C:79
 HESDfromKin.C:80
 HESDfromKin.C:81
 HESDfromKin.C:82
 HESDfromKin.C:83
 HESDfromKin.C:84
 HESDfromKin.C:85
 HESDfromKin.C:86
 HESDfromKin.C:87
 HESDfromKin.C:88
 HESDfromKin.C:89
 HESDfromKin.C:90
 HESDfromKin.C:91
 HESDfromKin.C:92
 HESDfromKin.C:93
 HESDfromKin.C:94
 HESDfromKin.C:95
 HESDfromKin.C:96
 HESDfromKin.C:97
 HESDfromKin.C:98
 HESDfromKin.C:99
 HESDfromKin.C:100
 HESDfromKin.C:101
 HESDfromKin.C:102
 HESDfromKin.C:103
 HESDfromKin.C:104
 HESDfromKin.C:105
 HESDfromKin.C:106
 HESDfromKin.C:107
 HESDfromKin.C:108
 HESDfromKin.C:109
 HESDfromKin.C:110
 HESDfromKin.C:111
 HESDfromKin.C:112
 HESDfromKin.C:113
 HESDfromKin.C:114
 HESDfromKin.C:115
 HESDfromKin.C:116
 HESDfromKin.C:117
 HESDfromKin.C:118
 HESDfromKin.C:119
 HESDfromKin.C:120
 HESDfromKin.C:121
 HESDfromKin.C:122
 HESDfromKin.C:123
 HESDfromKin.C:124
 HESDfromKin.C:125
 HESDfromKin.C:126
 HESDfromKin.C:127
 HESDfromKin.C:128
 HESDfromKin.C:129
 HESDfromKin.C:130
 HESDfromKin.C:131
 HESDfromKin.C:132
 HESDfromKin.C:133
 HESDfromKin.C:134
 HESDfromKin.C:135
 HESDfromKin.C:136
 HESDfromKin.C:137
 HESDfromKin.C:138