ROOT logo
TH1F *hpt;
void FASTSmear(char *finname="galice_upsi.root",char *foutname="FAST_upsi.root",Float_t bkg=0,Float_t nevmax=100000000){
 
  printf ("processing file %s\n",finname); 
  if (gClassTable->GetID("AliRun") < 0) {
    gROOT->LoadMacro("loadlibs.C");
    loadlibs();
  }    
  // create the structure to hold the variables for the branch
  
  struct fast_t {
    Float_t accp;
    Float_t effp;
    Float_t trigeffp;
    Float_t pgenp;
    Float_t psmearp; 
    Float_t thetagenp;
    Float_t thetasmearp;
    Float_t phigenp;
    Float_t phismearp;
    Float_t accm;
    Float_t effm;
    Float_t trigeffm;
    Float_t pgenm;
    Float_t psmearm; 
    Float_t thetagenm;
    Float_t thetasmearm;
    Float_t phigenm;
    Float_t phismearm;
  };

  fast_t fast;

  // open the input file
  TFile *fin=new TFile(finname);
  gAlice=(AliRun*)fin->Get("gAlice");

  // create the output file, its tree and branch 
  TFile *fout = new TFile(foutname,"RECREATE");
  TTree *FASTtrack = new TTree("FASTtrack","mu tracks smeared with fastsim");
  FASTtrack->Branch("fast",&fast.accp,"accp:effp:trigeffp:pgenp:psmearp:thetagenp:thetasmearp:phigenp:phismearp:accm:effm:trigeffm:pgenm:psmearm:thetagenm:thetasmearm:phigenm:phismearm");

  // create the fast tracker

  AliFastMuonTriggerEff *trigeff = new AliFastMuonTriggerEff();
  AliFastMuonTrackingAcc *acc = new AliFastMuonTrackingAcc();
  AliFastMuonTrackingEff *eff = new AliFastMuonTrackingEff();
  AliFastMuonTrackingRes *res = new AliFastMuonTrackingRes();
  acc->SetBackground(bkg);
  eff->SetBackground(bkg);
  res->SetBackground(bkg);  
  acc->Init(); 
  eff->Init(); 
  res->Init(); 
  AliFastDetector* tracker = new AliFastDetector("Tracker", "Muon Tracker");
  tracker->AddResponse(trigeff);
  tracker->AddResponse(acc);
  tracker->AddResponse(eff);
  tracker->AddResponse(res);
  tracker->Init();
  tracker->Dump();

  // loop over the events

  Int_t nev=AliRunLoader::GetNumberOfEvents();
  TParticle *mup, *mum;

  if (nev>nevmax) nev = nevmax;
  
  table = AliMUONFastTracking::Instance(); 
  printf ("background set to %g \n",table->GetBackground()); 
  for (Int_t iev=0; iev<nev; iev++) {
    Int_t npart=gAlice->GetEvent(iev);
    //    npart=500;
    for (Int_t ipart=0; ipart<npart; ipart+=3) {
      if (!(ipart%30)) printf ("Event #%d upsilon #%d \n",iev,ipart/3);
      part = gAlice->Particle(ipart+1);
      if (part->GetPdgCode() == 13) mum = part;
      else if (part->GetPdgCode() == -13) mup = part;
      part = gAlice->Particle(ipart+2);
      if (part->GetPdgCode() == 13) mum = part;
      else if (part->GetPdgCode() == -13) mup = part;

      // the mu+
      printf ("background set to %g \n",table->GetBackground()); 
      res->SetCharge(1);
      eff->SetCharge(1);
      acc->SetCharge(1);
      fast.pgenp     = mup->P();
      Double_t ptp    = fast.pgenp * TMath::Sin(mup->Theta());
      fast.thetagenp = 180./TMath::Pi() * mup->Theta();
      fast.phigenp   = 180./TMath::Pi() * mup->Phi();
      if (fast.phigenp>180) fast.phigenp -=360;
      res->Evaluate(fast.pgenp,fast.thetagenp,fast.phigenp,
		    fast.psmearp,fast.thetasmearp,fast.phismearp);
      fast.effp = eff->Evaluate(ptp,fast.thetagenp,fast.phigenp);
      fast.accp = acc->Evaluate(ptp,fast.thetagenp,fast.phigenp);
      fast.trigeffp = trigeff->Evaluate(1,ptp,fast.thetagenp,fast.phigenp);

      // the mu- 
      res->SetCharge(-1);
      acc->SetCharge(-1);
      eff->SetCharge(-1);
      fast.pgenm     = mum->P();
      Double_t ptm    = fast.pgenm * TMath::Sin(mum->Theta());
      fast.thetagenm = 180./TMath::Pi() * mum->Theta();
      fast.phigenm   = 180./TMath::Pi() * mum->Phi();
      if (fast.phigenm>180) fast.phigenm -=360;
      res->Evaluate(fast.pgenm,fast.thetagenm,fast.phigenm,
		    fast.psmearm,fast.thetasmearm,fast.phismearm);
      fast.effm = eff->Evaluate(ptm,fast.thetagenm,fast.phigenm);
      fast.accm = acc->Evaluate(ptm,fast.thetagenm,fast.phigenm);
      fast.trigeffm = trigeff->Evaluate(-1,ptm,fast.thetagenm,fast.phigenm);

      // fill the tree
      FASTtrack->Fill();
    }
  } 
  fout->Write();
  fout->Close();
  fin->Close();
}



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