ROOT logo
// Example: generation of kinematics tree with selected properties.
// Below we select events containing the decays D* -> D0 pi, D0 -> K- pi+
// inside the barrel part of the ALICE detector (45 < theta < 135)

#if !defined(__CINT__) || defined(__MAKECINT__)
#include <Riostream.h>
#include <TH1F.h>
#include <TStopwatch.h>
#include <TDatime.h>
#include <TRandom.h>
#include <TDatabasePDG.h>
#include <TParticle.h>
#include <TArrayI.h>
#include <TTree.h>

#include "AliGenerator.h"
#include "AliPDG.h"
#include "AliRunLoader.h"
#include "AliRun.h"
#include "AliStack.h"
#include "AliHeader.h"
#include "PYTHIA6/AliGenPythia.h"
#include "PYTHIA6/AliPythia.h"
#endif

Float_t EtaToTheta(Float_t arg);
void GetFinalDecayProducts(Int_t ind, AliStack & stack , TArrayI & ar);

void fastGen(Int_t nev = 1, const char* filename = "galice.root")
{
  AliPDG::AddParticlesToPdgDataBase();
  TDatabasePDG::Instance();
 


  // Run loader
  AliRunLoader* rl = AliRunLoader::Open("galice.root","FASTRUN","recreate");
  
  rl->SetCompressionLevel(2);
  rl->SetNumberOfEventsPerFile(nev);
  rl->LoadKinematics("RECREATE");
  rl->MakeTree("E");
  gAlice->SetRunLoader(rl);
  
  //  Create stack
  rl->MakeStack();
  AliStack* stack = rl->Stack();
  
  //  Header
  AliHeader* header = rl->GetHeader();
  
  //  Create and Initialize Generator
 
  // Example of charm generation taken from Config_PythiaHeavyFlavours.C
  AliGenPythia *gener = new AliGenPythia(-1);
  gener->SetEnergyCMS(14000.);
  gener->SetMomentumRange(0,999999);
  gener->SetPhiRange(0., 360.);
  gener->SetThetaRange(0.,180.);
  //  gener->SetProcess(kPyCharmppMNR); // Correct Pt distribution, wrong mult
  gener->SetProcess(kPyMb); // Correct multiplicity, wrong Pt
  gener->SetStrucFunc(kCTEQ4L);
  gener->SetPtHard(2.1,-1.0);
  gener->SetFeedDownHigherFamily(kFALSE);
  gener->SetStack(stack);
  gener->Init();

  // Go to galice.root
  rl->CdGAFile();

  // Forbid some decays. Do it after gener->Init(0, because
  // the initialization of the generator includes reading of the decay table.

  AliPythia * py= AliPythia::Instance();
  py->SetMDME(737,1,0); //forbid D*+->D+ + pi0
  py->SetMDME(738,1,0);//forbid D*+->D+ + gamma

  // Forbid all D0 decays except D0->K- pi+
  for(Int_t d=747; d<=762; d++){ 
    py->SetMDME(d,1,0);
  }
  // decay 763 is D0->K- pi+
  for(Int_t d=764; d<=807; d++){ 
    py->SetMDME(d,1,0);
  }

  //
  //                        Event Loop
  //
  
  TStopwatch timer;
  timer.Start();
  for (Int_t iev = 0; iev < nev; iev++) {
    
    cout <<"Event number "<< iev << endl;
    
    //  Initialize event
    header->Reset(0,iev);
    rl->SetEventNumber(iev);
    stack->Reset();
    rl->MakeTree("K");
    
    //  Generate event
    Int_t nprim = 0;
    Int_t ntrial = 0;
    Int_t ndstar = 0;
   
   

    //-------------------------------------------------------------------------------------

    while(!ndstar) {
      // Selection of events with D*
      stack->Reset();
      stack->ConnectTree(rl->TreeK());
      gener->Generate();
      ntrial++;
      nprim = stack->GetNprimary();
      
      for(Int_t ipart =0; ipart < nprim; ipart++){
        TParticle * part = stack->Particle(ipart);
        if(part)    {
          
          if (TMath::Abs(part->GetPdgCode())== 413) {

	    TArrayI daughtersId;

	    GetFinalDecayProducts(ipart,*stack,daughtersId);

	    Bool_t kineOK = kTRUE;

	    Double_t thetaMin = TMath::Pi()/4;
	    Double_t thetaMax = 3*TMath::Pi()/4;

	    for (Int_t id=1; id<=daughtersId[0]; id++) {
	      TParticle * daughter = stack->Particle(daughtersId[id]);
	      if (!daughter) {
		kineOK = kFALSE;
		break;
	      }

	      Double_t theta = daughter->Theta();
	      if (theta<thetaMin || theta>thetaMax) {
		kineOK = kFALSE;
		break;
	      }
	    }

	    if (!kineOK) continue;

            part->Print();
            ndstar++;     
	    
          }
        }
      }   
    }   
      
    cout << "Number of particles " << nprim << endl;
    cout << "Number of trials " << ntrial << endl;
    
    //  Finish event
    header->SetNprimary(stack->GetNprimary());
    header->SetNtrack(stack->GetNtrack());  
    
    //      I/O
    stack->FinishEvent();
    header->SetStack(stack);
    rl->TreeE()->Fill();
    rl->WriteKinematics("OVERWRITE");
    
  } // event loop
  timer.Stop();
  timer.Print();
  
  //                         Termination
  //  Generator
  gener->FinishRun();
  //  Write file
  rl->WriteHeader("OVERWRITE");
  gener->Write();
  rl->Write();
}



Float_t EtaToTheta(Float_t arg){
  return (180./TMath::Pi())*2.*atan(exp(-arg));
}


void GetFinalDecayProducts(Int_t ind, AliStack & stack , TArrayI & ar){

  // Recursive algorithm to get the final decay products of a particle
  //
  // ind is the index of the particle in the AliStack
  // stack is the particle stack from the generator
  // ar contains the indexes of the final decay products
  // ar[0] is the number of final decay products

  if (ind<0 || ind>stack.GetNtrack()) {
    cerr << "Invalid index of the particle " << ind << endl;
    return;
  } 
  if (ar.GetSize()==0) {
    ar.Set(10);
    ar[0] = 0;
  }

  TParticle * part = stack.Particle(ind);

  Int_t iFirstDaughter = part->GetFirstDaughter();
  if( iFirstDaughter<0) {
    // This particle is a final decay product, add its index to the array
    ar[0]++;
    if (ar.GetSize() <= ar[0]) ar.Set(ar.GetSize()+10); // resize if needed
    ar[ar[0]] = ind;
    return;
  } 

  Int_t iLastDaughter = part->GetLastDaughter();

  for (Int_t id=iFirstDaughter; id<=iLastDaughter;id++) {
    // Now search for final decay products of the daughters
    GetFinalDecayProducts(id,stack,ar);
  }
}
 fastGen.C:1
 fastGen.C:2
 fastGen.C:3
 fastGen.C:4
 fastGen.C:5
 fastGen.C:6
 fastGen.C:7
 fastGen.C:8
 fastGen.C:9
 fastGen.C:10
 fastGen.C:11
 fastGen.C:12
 fastGen.C:13
 fastGen.C:14
 fastGen.C:15
 fastGen.C:16
 fastGen.C:17
 fastGen.C:18
 fastGen.C:19
 fastGen.C:20
 fastGen.C:21
 fastGen.C:22
 fastGen.C:23
 fastGen.C:24
 fastGen.C:25
 fastGen.C:26
 fastGen.C:27
 fastGen.C:28
 fastGen.C:29
 fastGen.C:30
 fastGen.C:31
 fastGen.C:32
 fastGen.C:33
 fastGen.C:34
 fastGen.C:35
 fastGen.C:36
 fastGen.C:37
 fastGen.C:38
 fastGen.C:39
 fastGen.C:40
 fastGen.C:41
 fastGen.C:42
 fastGen.C:43
 fastGen.C:44
 fastGen.C:45
 fastGen.C:46
 fastGen.C:47
 fastGen.C:48
 fastGen.C:49
 fastGen.C:50
 fastGen.C:51
 fastGen.C:52
 fastGen.C:53
 fastGen.C:54
 fastGen.C:55
 fastGen.C:56
 fastGen.C:57
 fastGen.C:58
 fastGen.C:59
 fastGen.C:60
 fastGen.C:61
 fastGen.C:62
 fastGen.C:63
 fastGen.C:64
 fastGen.C:65
 fastGen.C:66
 fastGen.C:67
 fastGen.C:68
 fastGen.C:69
 fastGen.C:70
 fastGen.C:71
 fastGen.C:72
 fastGen.C:73
 fastGen.C:74
 fastGen.C:75
 fastGen.C:76
 fastGen.C:77
 fastGen.C:78
 fastGen.C:79
 fastGen.C:80
 fastGen.C:81
 fastGen.C:82
 fastGen.C:83
 fastGen.C:84
 fastGen.C:85
 fastGen.C:86
 fastGen.C:87
 fastGen.C:88
 fastGen.C:89
 fastGen.C:90
 fastGen.C:91
 fastGen.C:92
 fastGen.C:93
 fastGen.C:94
 fastGen.C:95
 fastGen.C:96
 fastGen.C:97
 fastGen.C:98
 fastGen.C:99
 fastGen.C:100
 fastGen.C:101
 fastGen.C:102
 fastGen.C:103
 fastGen.C:104
 fastGen.C:105
 fastGen.C:106
 fastGen.C:107
 fastGen.C:108
 fastGen.C:109
 fastGen.C:110
 fastGen.C:111
 fastGen.C:112
 fastGen.C:113
 fastGen.C:114
 fastGen.C:115
 fastGen.C:116
 fastGen.C:117
 fastGen.C:118
 fastGen.C:119
 fastGen.C:120
 fastGen.C:121
 fastGen.C:122
 fastGen.C:123
 fastGen.C:124
 fastGen.C:125
 fastGen.C:126
 fastGen.C:127
 fastGen.C:128
 fastGen.C:129
 fastGen.C:130
 fastGen.C:131
 fastGen.C:132
 fastGen.C:133
 fastGen.C:134
 fastGen.C:135
 fastGen.C:136
 fastGen.C:137
 fastGen.C:138
 fastGen.C:139
 fastGen.C:140
 fastGen.C:141
 fastGen.C:142
 fastGen.C:143
 fastGen.C:144
 fastGen.C:145
 fastGen.C:146
 fastGen.C:147
 fastGen.C:148
 fastGen.C:149
 fastGen.C:150
 fastGen.C:151
 fastGen.C:152
 fastGen.C:153
 fastGen.C:154
 fastGen.C:155
 fastGen.C:156
 fastGen.C:157
 fastGen.C:158
 fastGen.C:159
 fastGen.C:160
 fastGen.C:161
 fastGen.C:162
 fastGen.C:163
 fastGen.C:164
 fastGen.C:165
 fastGen.C:166
 fastGen.C:167
 fastGen.C:168
 fastGen.C:169
 fastGen.C:170
 fastGen.C:171
 fastGen.C:172
 fastGen.C:173
 fastGen.C:174
 fastGen.C:175
 fastGen.C:176
 fastGen.C:177
 fastGen.C:178
 fastGen.C:179
 fastGen.C:180
 fastGen.C:181
 fastGen.C:182
 fastGen.C:183
 fastGen.C:184
 fastGen.C:185
 fastGen.C:186
 fastGen.C:187
 fastGen.C:188
 fastGen.C:189
 fastGen.C:190
 fastGen.C:191
 fastGen.C:192
 fastGen.C:193
 fastGen.C:194
 fastGen.C:195
 fastGen.C:196
 fastGen.C:197
 fastGen.C:198
 fastGen.C:199
 fastGen.C:200
 fastGen.C:201
 fastGen.C:202
 fastGen.C:203
 fastGen.C:204
 fastGen.C:205
 fastGen.C:206
 fastGen.C:207
 fastGen.C:208
 fastGen.C:209
 fastGen.C:210
 fastGen.C:211
 fastGen.C:212
 fastGen.C:213
 fastGen.C:214
 fastGen.C:215
 fastGen.C:216
 fastGen.C:217
 fastGen.C:218
 fastGen.C:219
 fastGen.C:220
 fastGen.C:221
 fastGen.C:222
 fastGen.C:223
 fastGen.C:224
 fastGen.C:225
 fastGen.C:226
 fastGen.C:227
 fastGen.C:228