StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
testCPVDecays.cc
1 // Program to create ROOT files for EvtGen validation plots.
2 // This looks at the 1st generation daughters and stores 4-momenta
3 // info into a ROOT file for further analysis.
4 // Useful for Pythia, Photos and Tauola decay tests.
5 
6 #include "EvtGen/EvtGen.hh"
7 
8 #include "EvtGenBase/EvtParticle.hh"
9 #include "EvtGenBase/EvtParticleFactory.hh"
10 #include "EvtGenBase/EvtPatches.hh"
11 #include "EvtGenBase/EvtPDL.hh"
12 #include "EvtGenBase/EvtRandom.hh"
13 #include "EvtGenBase/EvtSimpleRandomEngine.hh"
14 #include "EvtGenBase/EvtMTRandomEngine.hh"
15 #include "EvtGenBase/EvtHepMCEvent.hh"
16 #include "EvtGenBase/EvtSpinDensity.hh"
17 #include "EvtGenBase/EvtSpinType.hh"
18 #include "EvtGenBase/EvtComplex.hh"
19 #include "EvtGenBase/EvtCPUtil.hh"
20 
21 #include "EvtGenBase/EvtAbsRadCorr.hh"
22 #include "EvtGenBase/EvtDecayBase.hh"
23 
24 #include "HepMC/GenEvent.h"
25 #include "HepMC/GenVertex.h"
26 #include "HepMC/GenParticle.h"
27 #include "HepMC/SimpleVector.h"
28 
29 #ifdef EVTGEN_EXTERNAL
30 #include "EvtGenExternal/EvtExternalGenList.hh"
31 #endif
32 
33 #include "TFile.h"
34 #include "TTree.h"
35 #include "TH1.h"
36 #include "TF1.h"
37 #include "TStyle.h"
38 #include "TCanvas.h"
39 #include "TLine.h"
40 #include "TROOT.h"
41 
42 #include <iostream>
43 #include <string>
44 #include <vector>
45 #include <list>
46 
47 using std::cout;
48 using std::endl;
49 using std::string;
50 
51 // Flight-time histograms for B0, B0bar
52 TH1F* H_total = new TH1F("Total", "", 300, 0.0, 12.0);
53 TH1F* H_B0 = new TH1F("B0", "", 300, 0.0, 12.0);
54 TH1F* H_B0bar = new TH1F("B0bar", "", 300, 0.0, 12.0);
55 
56 int B0Id(511), B0barId(-511);
57 
58 void storeBFlightTimes(HepMC::GenEvent* theEvent);
59 bool checkSignal(std::vector<int>& daugIdVect);
60 double calcFlightTime(HepMC::FourVector& BDecayVtx, HepMC::FourVector& B4mtm);
61 double sineFitFun(double* x, double* p);
62 double timeFitFun(double* x, double* p);
63 
64 int main(int argc, char** argv) {
65 
66  string decayFileName("CPVDecayFiles/Bd_JpsiKSeeCPV.dec");
67  if (argc > 1) {decayFileName = argv[1];}
68  cout<<"Decay file name is "<<decayFileName<<endl;
69 
70  string rootFileName("rootFiles/CPVDecayTest.root");
71  if (argc > 2) {rootFileName = argv[2];}
72  cout<<"Root file name is "<<rootFileName<<endl;
73 
74  string parentName("Upsilon(4S)");
75  if (argc > 3) {parentName = argv[3];}
76  cout<<"Parent name is "<<parentName<<endl;
77 
78  int nEvents(10000);
79  if (argc > 4) {nEvents = atoi(argv[4]);}
80 
81  double sin2Beta = sin(0.775);
82  if (argc > 5) {sin2Beta = atof(argv[5]);}
83 
84  cout<<"Number of events is "<<nEvents<<endl;
85  cout<<"sin2Beta = "<<sin2Beta<<" (used to draw oscillation maxima lines)"<<endl;
86 
87  // Define the random number generator
88 
89  EvtRandomEngine* myRandomEngine = 0;
90 
91 #ifdef EVTGEN_CPP11
92  // Use the Mersenne-Twister generator (C++11 only)
93  myRandomEngine = new EvtMTRandomEngine();
94 #else
95  myRandomEngine = new EvtSimpleRandomEngine();
96 #endif
97 
98  EvtAbsRadCorr* radCorrEngine = 0;
99  std::list<EvtDecayBase*> extraModels;
100 
101 #ifdef EVTGEN_EXTERNAL
102  EvtExternalGenList genList;
103  radCorrEngine = genList.getPhotosModel();
104  extraModels = genList.getListOfModels();
105 #endif
106 
107  int mixingType = EvtCPUtil::Incoherent;
108 
109  // Initialize the generator - read in the decay table and particle properties.
110  EvtGen myGenerator("../DECAY_2010.DEC", "../evt.pdl", myRandomEngine,
111  radCorrEngine, &extraModels, mixingType);
112 
113  myGenerator.readUDecay(decayFileName.c_str());
114 
115  EvtId theId = EvtPDL::getId(parentName);
116  if (theId.getId() == -1 && theId.getAlias() == -1) {
117  cout<<"Error. Could not find valid EvtId for "<<parentName<<endl;
118  return -1;
119  }
120 
121  // Start all initial (parent) decays at the origin
122  EvtVector4R origin(0.0, 0.0, 0.0, 0.0);
123  EvtSpinDensity* spinDensity = 0;
124 
125  EvtSpinType::spintype baseSpin = EvtPDL::getSpinType(theId);
126 
127  if (baseSpin == EvtSpinType::VECTOR) {
128  cout<<"Setting spin density for vector particle "<<parentName<<endl;
129  spinDensity = new EvtSpinDensity();
130  spinDensity->setDiag(EvtSpinType::getSpinStates(EvtSpinType::VECTOR));
131  spinDensity->set(1, 1, EvtComplex(0.0, 0.0));
132  }
133 
134  // Loop to create nEvents
135 
136  int i;
137  for (i = 0; i < nEvents; i++) {
138 
139  if (i%1000 == 0) {cout<<"Event number = "<<i+1<<" out of "<<nEvents<<std::endl;}
140 
141  // Set up the parent particle
142  int PDGId = EvtPDL::getStdHep(theId);
143  EvtVector4R pInit(EvtPDL::getMass(theId), 0.0, 0.0, 0.0);
144 
145  // Generate a new HepMC event with the decays. We own this pointer, so we must
146  // delete it after using the HepMC event information.
147  EvtHepMCEvent* theEvent = myGenerator.generateDecay(PDGId, pInit, origin, spinDensity);
148 
149  // Retrieve the HepMC event information
150  HepMC::GenEvent* hepMCEvent = theEvent->getEvent();
151  //hepMCEvent->print();
152 
153  // Fill the B0/B0bar flight time histograms
154  storeBFlightTimes(hepMCEvent);
155 
156  // Cleanup the event to avoid memory leaks
157  delete theEvent;
158 
159  }
160 
161  H_total->Sumw2();
162  H_B0->Sumw2();
163  H_B0bar->Sumw2();
164 
165  TH1F* H_Diff = dynamic_cast<TH1F*>(H_B0->Clone("H_Diff"));
166  H_Diff->Add(H_B0bar, -1.0);
167 
168  TH1F* H_DiffSum = dynamic_cast<TH1F*>(H_Diff->Clone("H_DiffSum"));
169  H_DiffSum->Divide(H_total);
170 
171  TF1* sineFit = new TF1("sineFit", sineFitFun, 0.0, 12.0, 3);
172 
173  sineFit->SetParName(0, "N");
174  sineFit->SetParName(1, "a");
175  sineFit->SetParName(2, "#phi");
176  sineFit->SetParameter(0, 0.5);
177  sineFit->SetParameter(1, 0.7);
178  sineFit->SetParameter(2, -0.7);
179  H_DiffSum->Fit(sineFit );
180 
181 
182  TF1* timeFit = new TF1("timeFit", timeFitFun, 0.0, 12.0, 2);
183  timeFit->SetParName(0, "N");
184  timeFit->SetParName(1, "#Gamma");
185  timeFit->SetParameter(0, 500);
186  timeFit->SetParameter(1, 0.6);
187  timeFit->SetParLimits(1, 0.0, 1.0);
188  H_total->Fit(timeFit);
189 
190  gROOT->SetStyle("Plain");
191  gStyle->SetOptFit(1111);
192  TCanvas* theCanvas = new TCanvas("theCanvas", "", 900, 700);
193  theCanvas->UseCurrentStyle();
194 
195  H_DiffSum->SetXTitle("t (ps)");
196  H_DiffSum->Draw();
197 
198  // Plot +- sin(2beta) lines
199  TLine line1(0.0, sin2Beta, 12.0, sin2Beta); line1.Draw();
200  TLine line2(0.0, -sin2Beta, 12.0, -sin2Beta); line2.Draw();
201  theCanvas->Print("BCPVSinFit.gif");
202 
203  H_total->SetXTitle("t (ps)");
204  H_total->Draw();
205  theCanvas->Print("BTimeFit.gif");
206 
207  TFile* theFile = new TFile(rootFileName.c_str(), "recreate");
208  theFile->cd();
209  H_B0->Write();
210  H_B0bar->Write();
211  H_total->Write();
212  H_Diff->Write();
213  H_DiffSum->Write();
214  theFile->Close();
215 
216  // Cleanup
217  delete theCanvas;
218  delete spinDensity;
219  delete myRandomEngine;
220 
221  cout<<"Done."<<endl;
222 
223  return 0;
224 
225 }
226 
227 void storeBFlightTimes(HepMC::GenEvent* theEvent) {
228 
229  std::list<HepMC::GenVertex*> allVertices;
231 
232  // Loop over vertices in the event
233  for (vertexIter = theEvent->vertices_begin();
234  vertexIter != theEvent->vertices_end(); ++vertexIter) {
235 
236  // Get the vertex
237  HepMC::GenVertex* theVertex = *vertexIter;
238 
239  if (theVertex == 0) {continue;}
240 
241  // Check to see if the incoming particle is a B candidate.
242  // If so, also look at the outgoing particles to see if we have a signal decay.
243  // For these, get the B decay vertex position and the B 4-momentum to calculate
244  // the B lifetime.
245 
246  bool gotB0(false), gotB0bar(false);
247  HepMC::FourVector B4mtm;
248 
250  for (inIter = theVertex->particles_in_const_begin();
251  inIter != theVertex->particles_in_const_end(); ++inIter) {
252 
253  HepMC::GenParticle* inParticle = *inIter;
254 
255  if (inParticle == 0) {continue;}
256 
257  int inPDGId = inParticle->pdg_id();
258  if (inPDGId == B0Id) {
259  gotB0 = true;
260  } else if (inPDGId == B0barId) {
261  gotB0bar = true;
262  }
263 
264  if (gotB0 == true || gotB0bar == true) {
265  B4mtm = inParticle->momentum();
266  }
267 
268  } // Loop over ingoing vertex particles
269 
270  if (gotB0 == true || gotB0bar == true) {
271 
272  // Check outgoing particles
273  std::vector<int> daugIdVect;
275  for (outIter = theVertex->particles_out_const_begin();
276  outIter != theVertex->particles_out_const_end(); ++outIter) {
277 
278  HepMC::GenParticle* outParticle = *outIter;
279 
280  if (outParticle != 0) {
281  int outPDGId = outParticle->pdg_id();
282  daugIdVect.push_back(outPDGId);
283  }
284 
285  } // Loop over outgoing vertex particles
286 
287  // Check if we have the signal decay
288  bool gotSignal = checkSignal(daugIdVect);
289 
290  // Fill the flight time histograms for signal B decays
291  if (gotSignal == true) {
292 
293  HepMC::FourVector BDecayVtx = theVertex->position();
294  double flightTime = calcFlightTime(BDecayVtx, B4mtm);
295 
296  if (gotB0 == true) {
297 
298  H_B0->Fill(flightTime);
299  H_total->Fill(flightTime);
300 
301  } else {
302 
303  H_B0bar->Fill(flightTime);
304  H_total->Fill(flightTime);
305 
306  }
307 
308  } // Got signal B decay (for flight-time histograms)
309 
310  } // Got a B candidate
311 
312  } // Loop over event vertices
313 
314 }
315 
316 double calcFlightTime(HepMC::FourVector& BDecayVtx, HepMC::FourVector& B4mtm) {
317 
318  double flightTime(0.0);
319 
320  double distance = BDecayVtx.rho()*1e-3; // metres
321  double momentum = B4mtm.rho(); // GeV/c
322  double BMass = 5.2795; // GeV/c^2
323  double c0 = 299792458.0; // m/s
324 
325  if (momentum > 0.0) {
326 
327  flightTime = 1.0e12*distance*BMass/(momentum*c0); // picoseconds
328 
329  }
330 
331  return flightTime;
332 
333 }
334 
335 bool checkSignal(std::vector<int>& daugIdVect) {
336 
337  bool gotSignal(false);
338 
339  int nDaug = daugIdVect.size();
340 
341  // Check for J/psi Ks decays
342  if (nDaug == 2) {
343 
344  int daug1Id = daugIdVect[0];
345  int daug2Id = daugIdVect[1];
346  if ((daug1Id == 443 && daug2Id == 310) || (daug1Id == 310 && daug2Id == 443)) {
347  gotSignal = true;
348  }
349 
350  }
351 
352  return gotSignal;
353 
354 }
355 
356 double sineFitFun(double* x, double* p) {
357 
358  double t = x[0];
359  double N = p[0];
360  double a = p[1];
361  double phi = p[2];
362 
363  double funcVal = N*sin(a*t + phi);
364 
365  return funcVal;
366 
367 }
368 
369 double timeFitFun(double* x, double* p) {
370 
371  double t = x[0];
372  double N0 = p[0];
373  double gamma = p[1];
374  double funcVal = N0*(exp(-gamma*t));
375 
376  return funcVal;
377 
378 }
int pdg_id() const
particle ID
Definition: GenParticle.h:214
std::vector< HepMC::GenParticle * >::const_iterator particles_in_const_iterator
const iterator for incoming particles
Definition: GenVertex.h:152
particles_out_const_iterator particles_out_const_end() const
end iteration of outgoing particles
Definition: GenVertex.h:450
std::vector< HepMC::GenParticle * >::const_iterator particles_out_const_iterator
const iterator for outgoing particles
Definition: GenVertex.h:155
GenVertex contains information about decay vertices.
Definition: GenVertex.h:52
double rho() const
spatial vector component magnitude
The GenEvent class is the core of HepMC.
Definition: GenEvent.h:155
Definition: EvtId.hh:27
particles_in_const_iterator particles_in_const_end() const
end iteration of incoming particles
Definition: GenVertex.h:440
vertex_const_iterator vertices_end() const
end vertex iteration
Definition: GenEvent.h:381
Definition: EvtGen.hh:46
particles_out_const_iterator particles_out_const_begin() const
begin iteration of outgoing particles
Definition: GenVertex.h:445
const FourVector & position() const
vertex position and time
Definition: GenVertex.h:406
const FourVector & momentum() const
standard 4 momentum
Definition: GenParticle.h:211
FourVector is a simple representation of a physics 4 vector.
Definition: SimpleVector.h:42
vertex_const_iterator vertices_begin() const
begin vertex iteration
Definition: GenEvent.h:377
non-const vertex iterator
Definition: GenEvent.h:391
particles_in_const_iterator particles_in_const_begin() const
begin iteration of incoming particles
Definition: GenVertex.h:435
The GenParticle class contains information about generated particles.
Definition: GenParticle.h:60