1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 /*****************************************************************
17 AliFlowOnTheFlyEventGenerator
19 Create events with realstic spectra and v_n values
20 store as AliFlowEventSimple
22 usage: see $ALICE_ROOT/PWGCF/FLOW/macros/GenerateOnTheFlyEvents.C
24 origin: Redmer Alexander Bertens, rbertens@cern.ch, rbertens@nikhef.nl, rbertens@uu.nl
25 Carlos Eugenio Perez Lara
27 ******************************************************************/
29 #include "Riostream.h"
30 #include "TObjArray.h"
31 #include "TClonesArray.h"
35 #include "TParticle.h"
43 #include "TParameter.h"
45 #include "AliFlowTrackSimple.h"
46 #include "AliFlowEventSimple.h"
48 #include "TVirtualMCDecayer.h"
49 #include "AliFlowOnTheFlyEventGenerator.h"
54 static Int_t OnTheFlyEventGeneratorCounter = 0;
56 ClassImp(AliFlowOnTheFlyEventGenerator)
57 //_____________________________________________________________________________
58 AliFlowOnTheFlyEventGenerator::AliFlowOnTheFlyEventGenerator() : fGenerators(0), fEmbedMe(0), fFlowEvent(0), fDecayer(0), fAddV2Mothers(0), fAddV3Mothers(0), fAddV2Daughters(0), fAddV3Daughters(0), fPsi2(0), fPsi3(0), fPrecisionPhi(.001), fMaxNumberOfIterations(100), fQA(0), fFF(0) {/* constructor used for IO by root */}
59 //_____________________________________________________________________________
60 AliFlowOnTheFlyEventGenerator::AliFlowOnTheFlyEventGenerator(Bool_t qa, Bool_t ff, Int_t mult, TVirtualMCDecayer* decayer, Bool_t a, Bool_t b, Bool_t c, Bool_t d) : fGenerators(0), fEmbedMe(0), fFlowEvent(0), fDecayer(0), fAddV2Mothers(0), fAddV3Mothers(0), fAddV2Daughters(0), fAddV3Daughters(0), fPsi2(0), fPsi3(0), fPrecisionPhi(.001), fMaxNumberOfIterations(100), fQA(qa), fFF(ff) {
62 if(!fFlowEvent) fFlowEvent = new AliFlowEventSimple(mult, (AliFlowEventSimple::ConstructionMethod)0, 0x0, 0, TMath::TwoPi(), -.9, .9);
64 fGenerators = new TObjArray();
65 fGenerators->SetOwner(kTRUE);
68 fDecayer = decayer; // decayer: user has ownership of decayer (see dtor)
73 if(fAddV3Mothers||fAddV3Daughters) printf(" > WARNING - v3 is not impelmented yet ! \n <");
75 //_____________________________________________________________________________
76 AliFlowOnTheFlyEventGenerator::~AliFlowOnTheFlyEventGenerator()
79 if(fGenerators) delete fGenerators;
80 if(fFlowEvent) delete fFlowEvent;
82 //_____________________________________________________________________________
83 AliFlowOnTheFlyEventGenerator::NaiveFlowAndSpectrumGenerator* AliFlowOnTheFlyEventGenerator::Find(Short_t pdg, Bool_t make)
85 // return instance of generator class for given fPDG value
86 for(int i(0); i < fGenerators->GetEntriesFast(); i++) {
87 if(((NaiveFlowAndSpectrumGenerator*)fGenerators->At(i))->GetPDGCode()==pdg) {
88 return (NaiveFlowAndSpectrumGenerator*)(fGenerators->At(i));
92 AliFlowOnTheFlyEventGenerator::NaiveFlowAndSpectrumGenerator* _tmp = new AliFlowOnTheFlyEventGenerator::NaiveFlowAndSpectrumGenerator(pdg, fQA, fFF);
93 fGenerators->Add(_tmp);
98 //_____________________________________________________________________________
99 void AliFlowOnTheFlyEventGenerator::SetPtSpectrum(const char* func, Short_t species)
101 // set pt spectrum for species, override default
102 NaiveFlowAndSpectrumGenerator* gen = Find(species, kTRUE);
103 TF1* spectrum = gen->GetPtSpectrum();
106 name = spectrum->GetName();
109 else name = Form("pt_%i", (int)species);
110 gen->SetPtSpectrum(new TF1(name.Data(), func, 0., 20));
112 //_____________________________________________________________________________
113 void AliFlowOnTheFlyEventGenerator::SetPtDependentV2(const char* func, Short_t species)
115 // set pt dependent v2 for species, override default
116 NaiveFlowAndSpectrumGenerator* gen = Find(species, kTRUE);
117 TF1* v2 = gen->GetDifferentialV2();
119 if(v2) { // if v2 is already present, get its name and release its memory
120 name = v2->GetName();
123 else name = Form("v2_%i", (int)species);
124 gen->SetDifferentialV2(new TF1(name.Data(),func, 0., 20));
126 //_____________________________________________________________________________
127 void AliFlowOnTheFlyEventGenerator::SetPtDependentV3(const char* func, Short_t species)
129 // set pt dependent v2 for species, override default
130 NaiveFlowAndSpectrumGenerator* gen = Find(species, kTRUE);
131 TF1* v3 = gen->GetDifferentialV3();
133 if(v3) { // if v3 is already present, get its name and release its memory
134 name = v3->GetName();
137 else name = Form("v3_%i", (int)species);
138 gen->SetDifferentialV3(new TF1(name.Data(), func, 0., 20));
140 //_____________________________________________________________________________
141 AliFlowEventSimple* AliFlowOnTheFlyEventGenerator::GenerateOnTheFlyEvent(TClonesArray* event, Int_t nSpecies, Int_t species[], Int_t mult[], Int_t bg, Bool_t fluc)
143 // generate a flow event according to settings provided in GENERATOR INTERFACE
145 fPsi2 = gRandom->Uniform(-TMath::Pi(), TMath::Pi());
146 // calculate the total multiplicity for the event
147 Int_t multPiKPr[] = {(int)(.7*bg/2.), (int)(.2*bg/2.), (int)(.1*bg/2.)};
148 Int_t fPDGPiKPr[] = {211, 321, 2212};
149 Int_t totalRPMultiplicity(0), totalPOIMultiplicity(0);
150 for(int i(0); i < nSpecies; i++) totalPOIMultiplicity+=mult[i];
151 for(int i(0); i < 3; i++) totalRPMultiplicity+=multPiKPr[i];
152 Int_t totalMultiplicity(totalRPMultiplicity+totalPOIMultiplicity);
153 // generate the particles of interest. if requested, a vn boost is given to the primary particles
154 for(int i(0); i < nSpecies; i++) {
155 if(fluc) GenerateOnTheFlyTracks(gRandom->Uniform(mult[i]-TMath::Sqrt(mult[i]), mult[i]+TMath::Sqrt(mult[i])), species[i], event);
156 else GenerateOnTheFlyTracks(mult[i], species[i], event);
158 // generate reference particles. if requested, a vn boost is given to the primary particles
159 for(int i(0); i < 3; i++) {
161 GenerateOnTheFlyTracks(gRandom->Uniform(multPiKPr[i]-TMath::Sqrt(multPiKPr[i]), multPiKPr[i]+TMath::Sqrt(mult[i])), fPDGPiKPr[i], event);
162 GenerateOnTheFlyTracks(gRandom->Uniform(multPiKPr[i]-TMath::Sqrt(multPiKPr[i]), multPiKPr[i]+TMath::Sqrt(mult[i])), -1*fPDGPiKPr[i], event);
165 GenerateOnTheFlyTracks(multPiKPr[i], fPDGPiKPr[i], event);
166 GenerateOnTheFlyTracks(multPiKPr[i], -1*fPDGPiKPr[i], event);
169 // if requested, decay the event
170 if(fDecayer) DecayOnTheFlyTracks(event);
171 // if an embedded event is given to the generator at this point, embed it
173 event->AbsorbObjects(fEmbedMe); // event has ownership of embedded event
174 fEmbedMe = 0x0; // reset to aviod embeddding more than once
176 // if requested, the secondaries are given a vn boost
177 if(fAddV2Daughters) AddV2(event);
178 // convert the event to an aliflowsimple event
179 // the tagging (rp) is done in this step, all tracks are tagged as rp's.
180 // daughters are made orphans (i.e. primaries with secondaries are rejected from the sample)
181 // converting clears the event tclones array
182 return ConvertTClonesToFlowEvent(event, totalMultiplicity);
184 //_____________________________________________________________________________
185 AliFlowEventSimple* AliFlowOnTheFlyEventGenerator::ConvertTClonesToFlowEvent(TClonesArray* event, Int_t totalMultiplicity)
187 // convert the content of a tclones array to an aliflowevent simple for output and tag particles as poi's and rp's
188 // first step, see if there's already a flow event available in memory. if not create one
189 if(!fFlowEvent) fFlowEvent = new AliFlowEventSimple(totalMultiplicity, (AliFlowEventSimple::ConstructionMethod)0, 0x0, 0, TMath::TwoPi(), -.9, .9);
190 // if there is a flow event, clear the members without releasing the allocated memory
191 else fFlowEvent->ClearFast();
192 fFlowEvent->SetMCReactionPlaneAngle(fPsi2);
194 TParticle* pParticle = 0x0;
195 Int_t iSelParticlesRP(0);
196 for (Int_t i(0); i < event->GetEntries(); i++) { // begin the loop over the input tracks
197 pParticle = (TParticle*)event->At(i); // get the particle
198 if (!pParticle) continue; // skip if empty slot (no particle)
199 if (pParticle->GetNDaughters()!=0) continue; // see if the particle has daughters (if so, reject it)
200 AliFlowTrackSimple* pTrack = new AliFlowTrackSimple(pParticle); // allocate space for the new flow track
201 pTrack->SetWeight(pParticle->Pz()); // ugly hack: store pz here ...
202 pTrack->SetID(pParticle->GetPdgCode()); // set pid code as id
203 pTrack->SetForRPSelection(kTRUE); // tag ALL particles as RP's,
204 // ugly hack 2: set charge to -1 for primaries and 1 for secondaries
205 if(pParticle->GetFirstMother()==-1) pTrack->SetCharge(-1);
206 else pTrack->SetCharge(1);
208 fFlowEvent->AddTrack(pTrack);
210 fFlowEvent->SetNumberOfRPs(iSelParticlesRP);
211 // all trakcs have safely been copied so we can clear the event
213 OnTheFlyEventGeneratorCounter = 0;
216 //_____________________________________________________________________________
217 void AliFlowOnTheFlyEventGenerator::AddV2(TParticle* part, Double_t v2)
219 // afterburner, adds v2, uses Newton-Raphson iteration
220 Double_t phi = part->Phi();
225 for (Int_t i=0; i!=fMaxNumberOfIterations; ++i) {
226 phiprev=phi; //store last value for comparison
227 f = phi-phi0+v2*TMath::Sin(2.*(phi-fPsi2));
228 fp = 1.0+2.0*v2*TMath::Cos(2.*(phi-fPsi2)); //first derivative
230 if (TMath::AreEqualAbs(phiprev,phi,fPrecisionPhi)) break;
232 part->SetMomentum( part->Pt()*TMath::Cos(phi), part->Pt()*TMath::Sin(phi), part->Pz(), part->Energy() );
234 //_____________________________________________________________________________
235 void AliFlowOnTheFlyEventGenerator::AddV2(TClonesArray* event)
237 // afterburner, adds v2 for different species to all tracks in an event
238 // FIXME at the moment no distincition between mothers and daughters
239 TParticle *part = new TParticle();
240 for(Int_t nTrack=0; nTrack!=event->GetEntriesFast(); ++nTrack) {
241 part = (TParticle*) event->At(nTrack);
242 // for each track, call overloaded addv2 function with the correct generator
243 // create a generator in the case where the decayer has introduced a new particle species
244 AddV2(part, Find(part->GetPdgCode(), kTRUE)->GetV2(part->Pt()));
247 //_____________________________________________________________________________
248 void AliFlowOnTheFlyEventGenerator::GenerateOnTheFlyTracks(Int_t mult, Int_t pid, TClonesArray* event)
250 // generate tracks for the event. the tracks that are generated here are the mothers
251 Double_t ph, et, pt, mm, px, py, pz, ee;
252 // access the descired generator through the static find method. if the particle species is unknown, we add it to the generator
253 NaiveFlowAndSpectrumGenerator* generator = Find(pid, kTRUE);
254 Int_t _tempCounter = OnTheFlyEventGeneratorCounter;
255 for(int i=_tempCounter; i<=mult+_tempCounter; ++i) {
256 OnTheFlyEventGeneratorCounter++;
257 TParticle* part = (TParticle*)event->ConstructedAt(i);
258 part->SetStatusCode(1);
259 part->SetFirstMother(-1);
260 part->SetLastMother(-1);
261 part->SetFirstDaughter(-1);
262 part->SetLastDaughter(-1);
263 ph = gRandom->Uniform(0,TMath::TwoPi());
264 et = gRandom->Uniform(-1.0,+1.0);
265 pt = generator->GetPt();
266 part->SetPdgCode( pid );
267 mm = part->GetMass();
268 px = pt*TMath::Cos(ph);
269 py = pt*TMath::Sin(ph);
270 pz = pt*TMath::SinH(et);
271 ee = TMath::Sqrt( mm*mm + pt*pt+pz*pz );
272 part->SetMomentum( px, py, pz, ee );
273 // if requested add v2 to the sample of mothers
274 if(fAddV2Mothers) AddV2(part, generator->GetV2(pt));
277 //_____________________________________________________________________________
278 void AliFlowOnTheFlyEventGenerator::DecayOnTheFlyTracks(TClonesArray *event)
280 // decay the tracks using a decayer that is set in the GenerateEventsOnTheFly.C macro
281 if(!fDecayer) return; // shouldn't happen ...
284 TClonesArray* arr = new TClonesArray("TParticle",10);
287 Int_t nTracks = event->GetEntriesFast();
288 TParticle *part, *part1;
289 TClonesArray &in = *event;
290 for(Int_t loop=0; loop!=3; ++loop) {
292 for(Int_t nTrack=nStart; nTrack!=nTracks; ++nTrack) {
294 part = (TParticle*) event->At( nTrack );
296 kf = TMath::Abs(part->GetPdgCode());
297 if(kf==22 ) ForceGammaDecay(arr, part); // if gamma, we decay it by hand
298 else { // else we let pythia do it
299 TLorentzVector pmom( part->Px(), part->Py(), part->Pz(), part->Energy() );
300 fDecayer->Decay(kf,&pmom);
301 fDecayer->ImportParticles(arr);
303 Int_t nDaughters = arr->GetEntries();
304 if(nDaughters<2) continue;
305 for(Int_t nDaughter=1; nDaughter!=nDaughters; ++nDaughter) {
306 part1 = (TParticle*) (arr->At(nDaughter));
307 if(!part1) continue; // was part1, mistake?
308 new(in[nTracks+secondaries]) TParticle( part1->GetPdgCode(),
309 part1->GetStatusCode(),
310 part1->GetFirstMother(),
311 part1->GetSecondMother(),
312 part1->GetFirstDaughter(),
313 part1->GetLastDaughter(),
323 if(nDaughter==1) part->SetFirstDaughter(nTracks+secondaries);
324 else if ((nDaughters-1)==nDaughter) part->SetLastDaughter(nTracks+secondaries);
325 else part->SetDaughter(nDaughter,nTracks+secondaries);
329 nTracks = event->GetEntries();
334 //_____________________________________________________________________________
335 void AliFlowOnTheFlyEventGenerator::ForceGammaDecay(TClonesArray* arr, TParticle* part)
338 const Float_t mass1(gRandom->Uniform(0.0000001, .1)); // FIXME randomly sample photon 'mass'
339 const Float_t mass_d1(.000511); // daughter mass (electron)
340 const Float_t mass_d2(.000511);
341 Double_t p_dec1=sqrt((TMath::Abs(mass1*mass1-(mass_d1+mass_d2)*(mass_d1+mass_d2))*(mass1*mass1-(mass_d1-mass_d2)*(mass_d1-mass_d2))))/2/mass1; // energy
342 TLorentzVector *p_d1=new TLorentzVector(); // lorentz vectors for daughters
343 TLorentzVector *p_d2=new TLorentzVector();
344 Double_t yy_m=part->Y(); // pseudo rapidity of mother
345 Double_t pt_m=part->Pt(); // pt of mother
346 Double_t mt_m=sqrt(mass1*mass1+pt_m*pt_m); // transverse mass of mother
347 Double_t pz_m=mt_m*TMath::SinH(yy_m); // MAGIC
348 Double_t phi_m=gRandom->Rndm()*2*TMath::Pi();
349 Double_t px_m=pt_m*TMath::Cos(phi_m);
350 Double_t py_m=pt_m*TMath::Sin(phi_m);
351 Double_t costh_d=2*gRandom->Rndm()-1;
352 Double_t phi_d=gRandom->Rndm()*2*TMath::Pi();
353 Double_t pz_d=p_dec1*costh_d;
354 Double_t pt_d=p_dec1*sqrt(1-costh_d*costh_d);
355 Double_t px_d=pt_d*TMath::Cos(phi_d);
356 Double_t py_d=pt_d*TMath::Sin(phi_d);
357 p_d1->SetPxPyPzE(px_d,py_d,pz_d,sqrt(mass_d1*mass_d1+p_dec1*p_dec1));
358 p_d2->SetPxPyPzE(-px_d,-py_d,-pz_d,sqrt(mass_d2*mass_d2+p_dec1*p_dec1));
359 Double_t gamma_b=sqrt(mass1*mass1+pz_m*pz_m+pt_m*pt_m)/mass1;
360 Double_t bx=px_m/gamma_b/mass1;
361 Double_t by=py_m/gamma_b/mass1;
362 Double_t bz=pz_m/gamma_b/mass1;
363 p_d1->Boost(bx,by,bz);
364 p_d2->Boost(bx,by,bz);
365 TParticle* daughter_1 = (TParticle*)arr->ConstructedAt(0);
366 TParticle* daughter_2 = (TParticle*)arr->ConstructedAt(1);
367 daughter_1->SetStatusCode(1);
368 daughter_1->SetFirstMother(-1);
369 daughter_1->SetLastMother(-1);
370 daughter_1->SetFirstDaughter(-1);
371 daughter_1->SetLastDaughter(-1);
372 daughter_1->SetPdgCode(11);
373 TLorentzVector& ref_p_d1 = *p_d1;
374 daughter_1->SetMomentum(ref_p_d1);
375 daughter_2->SetStatusCode(1);
376 daughter_2->SetFirstMother(-1);
377 daughter_2->SetLastMother(-1);
378 daughter_2->SetFirstDaughter(-1);
379 daughter_2->SetLastDaughter(-1);
380 daughter_2->SetPdgCode(-11);
381 TLorentzVector& pp = *p_d2;
382 daughter_2->SetMomentum(pp);
384 //_____________________________________________________________________________
385 void AliFlowOnTheFlyEventGenerator::InitGenerators()
387 // initializes generators for a list of common particles
388 // new generators will be created 'on the fly' if they are encountered
389 // if a generator is requested for a non-existent pdg value, things will get unstable ...
390 Short_t PDGcode[] = {11, -11, // e-
418 2212, -2212, // proton
420 3312, -3312, // xi- (cascade)
425 for(int i(0); i < 47; i++) fGenerators->Add(new NaiveFlowAndSpectrumGenerator(PDGcode[i], fQA, fFF));
427 //_____________________________________________________________________________
428 void AliFlowOnTheFlyEventGenerator::PrintGenerators()
430 // make sure all went ok
432 printf(" > PANIC <\n");
435 printf(" > %i species available in generator\n", fGenerators->GetEntriesFast());
436 for(int i(0); i < fGenerators->GetEntriesFast(); i++)
437 printf(" - PDG: %i \n", ((NaiveFlowAndSpectrumGenerator*)fGenerators->At(i))->GetPDGCode());
438 printf(" more can be created 'on the fly', but beware of non-existent pdg values !");
440 //_____________________________________________________________________________
441 void AliFlowOnTheFlyEventGenerator::DoGeneratorQA(Bool_t v2, Bool_t v3)
443 // this function loops over all the species that are available in the fGenerators, which are
444 // a) created by the user (in the event generator) as mother particles or
445 // b) introduced by pythia (the decayer) as daughter particles
446 // c) made by the InitGenerators() function
447 // 'empty species' (i.e. species for which a generator exists but which were not actually sampled) are omitted
448 // the output histograms are written to a file named "GeneratorfQA.root"
450 printf(" > Request has been made for QA plots but QA histograms have been disabled !\n");
453 TFile *QAfile = new TFile("GeneratorfQA.root","RECREATE");
454 for(int i(0); i < fGenerators->GetEntriesFast(); i++) {
455 NaiveFlowAndSpectrumGenerator* gen = (NaiveFlowAndSpectrumGenerator*)fGenerators->At(i);
456 if(!gen) continue; // skip if we failed to get a generator
457 TH1F* QApt = (TH1F*)gen->GetQAType(0);
458 TH2F* QAv2 = (TH2F*)gen->GetQAType(1);
459 TH2F* QAv3 = (TH2F*)gen->GetQAType(2);
460 if((!QApt)||(!QAv2)||(!QAv3)) {
461 printf(" > couldn't read qa histogrmas for species %i <\n", gen->GetPDGCode());
464 if(QApt->GetEntries()==0&&QAv2->GetEntries()==0&&QAv3->GetEntries()==0) continue; // skip if no tracks of this species have been sampled
465 printf(" > saving QA histograms for sampled species %i <\n", gen->GetPDGCode());
466 QAfile->mkdir(Form("fPDG_%i", gen->GetPDGCode())); // create a directory in the output file
467 QAfile->cd(Form("fPDG_%i", gen->GetPDGCode())); // cd into this directory
468 if(!QApt->GetEntries()==0) { // only plot the pt fSpectrum if this guy was generated
469 if(!gen->GetPtSpectrum()->Integral(0,20)==0) {
470 QApt->Scale(gen->GetPtSpectrum()->Integral(0,20)/QApt->Integral(),"width");
474 gen->GetPtSpectrum()->SetNpx(10000); // otherwise tf1 plot is very ugly
475 gen->GetPtSpectrum()->Draw();
476 gen->GetPtSpectrum()->Write();
479 gen->GetDifferentialV2()->SetNpx(10000);
480 gen->GetDifferentialV2()->Draw();
481 gen->GetDifferentialV2()->Write();
486 gen->GetDifferentialV3()->Draw();
487 gen->GetDifferentialV3()->SetNpx(10000);
488 gen->GetDifferentialV3()->Write();
494 //_____________________________________________________________________________