CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TauValidation.cc
Go to the documentation of this file.
1 /*class TauValidation
2  *
3  * Class to fill dqm monitor elements from existing EDM file
4  *
5  */
7 
8 #include "CLHEP/Units/defs.h"
9 #include "CLHEP/Units/PhysicalConstants.h"
10 
12 
16 
17 using namespace edm;
18 
20  wmanager_(iPSet,consumesCollector())
21  ,hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection"))
22  ,tauEtCut(iPSet.getParameter<double>("tauEtCutForRtau"))
23  ,NJAKID(22)
24  ,zsbins(20)
25  ,zsmin(-0.5)
26  ,zsmax(0.5)
27 {
28  dbe = 0;
30 
31  hepmcCollectionToken_=consumes<HepMCProduct>(hepmcCollection_);
32 }
33 
35 
37 {
38  if(dbe){
40  dbe->setCurrentFolder("Generator/Tau");
41 
42  // Number of analyzed events
43  nTaus = dbe->book1D("nTaus", "n analyzed Taus", 1, 0., 1.);
44  nPrimeTaus = dbe->book1D("nPrimeTaus", "n analyzed prime Taus", 1, 0., 1.);
45 
46  //Kinematics
47  TauPt = dbe->book1D("TauPt","Tau pT", 100 ,0,100);
48  TauEta = dbe->book1D("TauEta","Tau eta", 100 ,-2.5,2.5);
49  TauPhi = dbe->book1D("TauPhi","Tau phi", 100 ,-3.14,3.14);
50  TauProngs = dbe->book1D("TauProngs","Tau n prongs", 7 ,0,7);
51  TauDecayChannels = dbe->book1D("TauDecayChannels","Tau decay channels", 13 ,0,13);
55  TauDecayChannels->setBinLabel(1+pi,"#pi^{#pm}");
56  TauDecayChannels->setBinLabel(1+rho,"#rho^{#pm}");
57  TauDecayChannels->setBinLabel(1+a1,"a_{1}^{#pm}");
58  TauDecayChannels->setBinLabel(1+pi1pi0,"#pi^{#pm}#pi^{0}");
59  TauDecayChannels->setBinLabel(1+pinpi0,"#pi^{#pm}n#pi^{0}");
60  TauDecayChannels->setBinLabel(1+tripi,"3#pi^{#pm}");
61  TauDecayChannels->setBinLabel(1+tripinpi0,"3#pi^{#pm}n#pi^{0}");
64  TauDecayChannels->setBinLabel(1+stable,"Stable");
65 
66  TauMothers = dbe->book1D("TauMothers","Tau mother particles", 10 ,0,10);
68  TauMothers->setBinLabel(1+B,"B Decays");
69  TauMothers->setBinLabel(1+D,"D Decays");
70  TauMothers->setBinLabel(1+gamma,"#gamma");
71  TauMothers->setBinLabel(1+Z,"Z");
72  TauMothers->setBinLabel(1+W,"W");
73  TauMothers->setBinLabel(1+HSM,"H_{SM}/h^{0}");
74  TauMothers->setBinLabel(1+H0,"H^{0}");
75  TauMothers->setBinLabel(1+A0,"A^{0}");
76  TauMothers->setBinLabel(1+Hpm,"H^{#pm}");
77 
78  DecayLength = dbe->book1D("DecayLength","#tau Decay Length", 100 ,0,20); DecayLength->setAxisTitle("L_{#tau} (mm)");
79  LifeTime = dbe->book1D("LifeTime","#tau LifeTime ", 100 ,0,1000E-15); LifeTime->setAxisTitle("#tau_{#tau}s)");
80 
81  TauRtauW = dbe->book1D("TauRtauW","W->Tau p(leading track)/E(visible tau)", 50 ,0,1); TauRtauW->setAxisTitle("rtau");
82  TauRtauHpm = dbe->book1D("TauRtauHpm","Hpm->Tau p(leading track)/E(visible tau)", 50 ,0,1); TauRtauHpm->setAxisTitle("rtau");
83 
84  TauSpinEffectsW_X = dbe->book1D("TauSpinEffectsWX","Pion energy in W rest frame", 50 ,0,1); TauSpinEffectsW_X->setAxisTitle("X");
85  TauSpinEffectsHpm_X = dbe->book1D("TauSpinEffectsHpmX","Pion energy in Hpm rest frame", 50 ,0,1); TauSpinEffectsHpm_X->setAxisTitle("X");
86 
87  TauSpinEffectsW_eX = dbe->book1D("TauSpinEffectsWeX","e energy in W rest frame", 50 ,0,1); TauSpinEffectsW_eX->setAxisTitle("X");
88  TauSpinEffectsHpm_eX = dbe->book1D("TauSpinEffectsHpmeX","e energy in Hpm rest frame", 50 ,0,1); TauSpinEffectsHpm_eX->setAxisTitle("X");
89 
90  TauSpinEffectsW_muX = dbe->book1D("TauSpinEffectsWmuX","mu energy in W rest frame", 50 ,0,1); TauSpinEffectsW_muX->setAxisTitle("X");
91  TauSpinEffectsHpm_muX = dbe->book1D("TauSpinEffectsHpmmuX","mu energy in Hpm rest frame", 50 ,0,1); TauSpinEffectsHpm_muX->setAxisTitle("X");
92 
93  TauSpinEffectsW_UpsilonRho = dbe->book1D("TauSpinEffectsWUpsilonRho","#Upsilon for #rho", 50 ,-1,1); TauSpinEffectsW_UpsilonRho->setAxisTitle("#Upsilon");
94  TauSpinEffectsHpm_UpsilonRho = dbe->book1D("TauSpinEffectsHpmUpsilonRho","#Upsilon for #rho", 50 ,-1,1); TauSpinEffectsHpm_UpsilonRho->setAxisTitle("#Upsilon");
95 
96  TauSpinEffectsW_UpsilonA1 = dbe->book1D("TauSpinEffectsWUpsilonA1","#Upsilon for a1", 50 ,-1,1); TauSpinEffectsW_UpsilonA1->setAxisTitle("#Upsilon");
97  TauSpinEffectsHpm_UpsilonA1 = dbe->book1D("TauSpinEffectsHpmUpsilonA1","#Upsilon for a1", 50 ,-1,1); TauSpinEffectsHpm_UpsilonA1->setAxisTitle("#Upsilon");
98 
99  TauSpinEffectsH_pipiAcoplanarity = dbe->book1D("TauSpinEffectsH_pipiAcoplanarity","H Acoplanarity for #pi^{-}#pi^{+}", 50 ,0,2*TMath::Pi()); TauSpinEffectsH_pipiAcoplanarity->setAxisTitle("Acoplanarity");
100 
101  TauSpinEffectsH_pipiAcollinearity = dbe->book1D("TauSpinEffectsH_pipiAcollinearity","H Acollinearity for #pi^{-}#pi^{+}", 50 ,0,TMath::Pi()); TauSpinEffectsH_pipiAcollinearity->setAxisTitle("Acollinearity");
102  TauSpinEffectsH_pipiAcollinearityzoom = dbe->book1D("TauSpinEffectsH_pipiAcollinearityzoom","H Acollinearity for #pi^{-}#pi^{+}", 50 ,3,TMath::Pi()); TauSpinEffectsH_pipiAcollinearityzoom->setAxisTitle("Acollinearity");
103 
104  TauSpinEffectsZ_MVis = dbe->book1D("TauSpinEffectsZMVis","Mass of pi+ pi-", 25 ,0,1.1); TauSpinEffectsZ_MVis->setAxisTitle("M_{#pi^{+}#pi^{-}}");
105  TauSpinEffectsH_MVis = dbe->book1D("TauSpinEffectsHMVis","Mass of pi+ pi-", 25 ,0,1.1); TauSpinEffectsZ_MVis->setAxisTitle("M_{#pi^{+}#pi^{-}}");
106 
107  TauSpinEffectsZ_Zs = dbe->book1D("TauSpinEffectsZZs","Z_{s}", zsbins ,zsmin,zsmax); TauSpinEffectsZ_Zs->setAxisTitle("Z_{s}");
108  TauSpinEffectsH_Zs = dbe->book1D("TauSpinEffectsHZs","Z_{s}", zsbins ,zsmin,zsmax); TauSpinEffectsZ_Zs->setAxisTitle("Z_{s}");
109 
110  TauSpinEffectsZ_X= dbe->book1D("TauSpinEffectsZX","X of #tau^{-}", 25 ,0,1.0); TauSpinEffectsZ_X->setAxisTitle("X");
111  TauSpinEffectsH_X= dbe->book1D("TauSpinEffectsH_X","X of #tau^{-}", 25 ,0,1.0); TauSpinEffectsH_X->setAxisTitle("X");
112 
113  TauSpinEffectsZ_Xf = dbe->book1D("TauSpinEffectsZXf","X of forward emitted #tau^{-}", 25 ,0,1.0); TauSpinEffectsZ_Xf->setAxisTitle("X_{f}");
114  TauSpinEffectsH_Xf = dbe->book1D("TauSpinEffectsHXf","X of forward emitted #tau^{-}", 25 ,0,1.0); TauSpinEffectsZ_Xf->setAxisTitle("X_{f}");
115 
116  TauSpinEffectsZ_Xb = dbe->book1D("TauSpinEffectsZXb","X of backward emitted #tau^{-}", 25 ,0,1.0); TauSpinEffectsZ_Xb->setAxisTitle("X_{b}");
117  TauSpinEffectsH_Xb = dbe->book1D("TauSpinEffectsHXb","X of backward emitted #tau^{-}", 25 ,0,1.0); TauSpinEffectsZ_Xb->setAxisTitle("X_{b}");
118 
119  TauSpinEffectsZ_eX = dbe->book1D("TauSpinEffectsZeX","e energy in Z rest frame", 50 ,0,1); TauSpinEffectsZ_eX->setAxisTitle("X");
120  TauSpinEffectsH_eX = dbe->book1D("TauSpinEffectsHeX","e energy in H rest frame", 50 ,0,1); TauSpinEffectsH_eX->setAxisTitle("X");
121 
122  TauSpinEffectsZ_muX = dbe->book1D("TauSpinEffectsZmuX","mu energy in z rest frame", 50 ,0,1); TauSpinEffectsZ_muX->setAxisTitle("X");
123  TauSpinEffectsH_muX = dbe->book1D("TauSpinEffectsHmuX","mu energy in H rest frame", 50 ,0,1); TauSpinEffectsH_muX->setAxisTitle("X");
124 
125  TauSpinEffectsH_rhorhoAcoplanarityminus = dbe->book1D("TauSpinEffectsH_rhorhoAcoplanarityminus","#phi^{*-} (acoplanarity) for Higgs #rightarrow #rho-#rho (y_{1}*y_{2}<0)", 32 ,0,2*TMath::Pi()); TauSpinEffectsH_rhorhoAcoplanarityminus->setAxisTitle("#phi^{*-} (Acoplanarity)");
126  TauSpinEffectsH_rhorhoAcoplanarityplus = dbe->book1D("TauSpinEffectsH_rhorhoAcoplanarityplus","#phi^{*+} (acoplanarity) for Higgs #rightarrow #rho-#rho (y_{1}*y_{2}>0)", 32 ,0,2*TMath::Pi()); TauSpinEffectsH_rhorhoAcoplanarityplus->setAxisTitle("#phi^{*+} (Acoplanarity)");
127 
128  TauFSRPhotonsN=dbe->book1D("TauFSRPhotonsN","FSR Photons radiating from/with tau (Gauge Boson)", 5 ,-0.5,4.5);
129  TauFSRPhotonsN->setAxisTitle("N FSR Photons radiating from/with tau");
130  TauFSRPhotonsPt=dbe->book1D("TauFSRPhotonsPt","Pt of FSR Photons radiating from/with tau (Gauge Boson)", 100 ,0,100);
131  TauFSRPhotonsPt->setAxisTitle("P_{t} of FSR Photons radiating from/with tau [per tau]");
132  TauFSRPhotonsPtSum=dbe->book1D("TauFSRPhotonsPtSum","Pt of FSR Photons radiating from/with tau (Gauge Boson)", 100 ,0,100);
133  TauFSRPhotonsPtSum->setAxisTitle("P_{t} of FSR Photons radiating from/with tau [per tau]");
134 
135  TauBremPhotonsN=dbe->book1D("TauBremPhotonsN","Brem. Photons radiating in tau decay", 5 ,-0.5,4.5);
136  TauBremPhotonsN->setAxisTitle("N FSR Photons radiating from/with tau");
137  TauBremPhotonsPt=dbe->book1D("TauBremPhotonsPt","Sum Brem Pt ", 100 ,0,100);
138  TauFSRPhotonsPt->setAxisTitle("P_{t} of Brem. Photons radiating in tau decay");
139  TauBremPhotonsPtSum =dbe->book1D("TauBremPhotonsPtSum","Sum of Brem Pt ", 100 ,0,100);
140  TauFSRPhotonsPtSum->setAxisTitle("Sum P_{t} of Brem. Photons radiating in tau decay");
141 
142  JAKID =dbe->book1D("JAKID","JAK ID",NJAKID+1,-0.5,NJAKID+0.5);
143  for(unsigned int i=0; i<NJAKID+1;i++){
144  JAKInvMass.push_back(std::vector<MonitorElement *>());
145  TString tmp="JAKID";
146  tmp+=i;
147  JAKInvMass.at(i).push_back(dbe->book1D("M"+tmp,"M_{"+tmp+"} (GeV)", 80 ,0,2.0));
148  if(i==TauDecay::JAK_A1_3PI ||
151  JAKInvMass.at(i).push_back(dbe->book1D("M13"+tmp,"M_{13,"+tmp+"} (GeV)", 80 ,0,2.0));
152  JAKInvMass.at(i).push_back(dbe->book1D("M23"+tmp,"M_{23,"+tmp+"} (GeV)", 80 ,0,2.0));
153  JAKInvMass.at(i).push_back(dbe->book1D("M12"+tmp,"M_{12,"+tmp+"} (GeV)", 80 ,0,2.0));
154  }
155  }
156  }
157 
158  return;
159 }
160 
162  return;
163 }
164 
165 void TauValidation::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup)
166 {
168  iSetup.getData( fPDGTable );
169  return;
170 }
171 void TauValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
173 {
176  iEvent.getByToken(hepmcCollectionToken_, evt);
177 
178  //Get EVENT
179  HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
180 
181  double weight = wmanager_.weight(iEvent);
182 
184  /*
185  edm::Handle<double> WT;
186  iEvent.getByLabel(edm::InputTag("TauSpinnerGen","TauSpinnerWT"),WT);
187  weight = 1.0;
188  if(*(WT.product())>1e-3 && *(WT.product())<=10.0) weight=(*(WT.product()));
189  else {weight=1.0;}
190  */
192 
193  // find taus
194  for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); iter++) {
195  if(abs((*iter)->pdg_id())==PdtPdgMini::Z0 || abs((*iter)->pdg_id())==PdtPdgMini::Higgs0){
196  spinEffectsZH(*iter,weight);
197  }
198  if(abs((*iter)->pdg_id())==15){
199  if(isLastTauinChain(*iter)){
200  nTaus->Fill(0.5,weight);
201  int mother = tauMother(*iter,weight);
202  int decaychannel = tauDecayChannel(*iter,weight);
203  tauProngs(*iter, weight);
204  if(mother>-1){ // exclude B, D and other non-signal decay modes
205  nPrimeTaus->Fill(0.5,weight);
206  TauPt->Fill((*iter)->momentum().perp(),weight);
207  TauEta->Fill((*iter)->momentum().eta(),weight);
208  TauPhi->Fill((*iter)->momentum().phi(),weight);
209  rtau(*iter,mother,decaychannel,weight);
210  photons(*iter,weight);
211  }
213  //Adding JAKID and Mass information
214  //
215  TauDecay_CMSSW TD;
216  unsigned int jak_id, TauBitMask;
217  if(TD.AnalyzeTau((*iter),jak_id,TauBitMask,false,false)){
218  JAKID->Fill(jak_id,weight);
219  if(jak_id<=NJAKID){
220  int tcharge=(*iter)->pdg_id()/abs((*iter)->pdg_id());
221  std::vector<HepMC::GenParticle*> part=TD.Get_TauDecayProducts();
222  spinEffectsWHpm(*iter,mother,jak_id,part,weight);
223  TLorentzVector LVQ(0,0,0,0);
224  TLorentzVector LVS12(0,0,0,0);
225  TLorentzVector LVS13(0,0,0,0);
226  TLorentzVector LVS23(0,0,0,0);
227  bool haspart1=false;
228  TVector3 PV,SV;
229  bool hasDL(false);
230  for(unsigned int i=0;i<part.size();i++){
231  if(abs(part.at(i)->pdg_id())!=PdtPdgMini::nu_tau && TD.isTauFinalStateParticle(part.at(i)->pdg_id()) && !hasDL){
232  PV=TVector3((*iter)->production_vertex()->point3d().x(),(*iter)->production_vertex()->point3d().y(),(*iter)->production_vertex()->point3d().z());
233  SV=TVector3(part.at(i)->production_vertex()->point3d().x(),part.at(i)->production_vertex()->point3d().y(),part.at(i)->production_vertex()->point3d().z());
234  TVector3 DL=SV-PV;
235  DecayLength->Fill(DL.Mag(),weight);
236  double c(2.99792458E8),Ltau(DL.Mag()/1000),beta((*iter)->momentum().rho()/(*iter)->momentum().m());
237  LifeTime->Fill( Ltau/(c*beta),weight);
238  }
239 
240  if(TD.isTauFinalStateParticle(part.at(i)->pdg_id()) &&
241  abs(part.at(i)->pdg_id())!=PdtPdgMini::nu_e &&
242  abs(part.at(i)->pdg_id())!=PdtPdgMini::nu_mu &&
243  abs(part.at(i)->pdg_id())!=PdtPdgMini::nu_tau ){
244  TLorentzVector LV(part.at(i)->momentum().px(),part.at(i)->momentum().py(),part.at(i)->momentum().pz(),part.at(i)->momentum().e());
245  LVQ+=LV;
246  if(jak_id==TauDecay::JAK_A1_3PI ||
247  jak_id==TauDecay::JAK_KPIK ||
248  jak_id==TauDecay::JAK_KPIPI
249  ){
250  if((tcharge==part.at(i)->pdg_id()/abs(part.at(i)->pdg_id()) && TD.nProng(TauBitMask)==3) || (jak_id==TauDecay::JAK_A1_3PI && TD.nProng(TauBitMask)==1 && abs(part.at(i)->pdg_id())==PdtPdgMini::pi_plus) ){
251  LVS13+=LV;
252  LVS23+=LV;
253  }
254  else{
255  LVS12+=LV;
256  if(!haspart1 && ((jak_id==TauDecay::JAK_A1_3PI) || (jak_id!=TauDecay::JAK_A1_3PI && abs(part.at(i)->pdg_id())==PdtPdgMini::K_plus) )){
257  LVS13+=LV;
258  haspart1=true;
259  }
260  else{
261  LVS23+=LV;
262  }
263  }
264  }
265  }
266  }
267  part.clear();
268  JAKInvMass.at(jak_id).at(0)->Fill(LVQ.M(),weight);
269  if(jak_id==TauDecay::JAK_A1_3PI ||
270  jak_id==TauDecay::JAK_KPIK ||
271  jak_id==TauDecay::JAK_KPIPI
272  ){
273  JAKInvMass.at(jak_id).at(1)->Fill(LVS13.M(),weight);
274  JAKInvMass.at(jak_id).at(2)->Fill(LVS23.M(),weight);
275  JAKInvMass.at(jak_id).at(3)->Fill(LVS12.M(),weight);
276  }
277  }
278  }
279  else{
280  JAKID->Fill(jak_id,weight);
281  }
282  }
283  }
284  }
285  delete myGenEvent;
286 }//analyze
287 
289  if ( tau->production_vertex() ) {
290  HepMC::GenVertex::particle_iterator mother;
291  for (mother = tau->production_vertex()->particles_begin(HepMC::parents); mother!= tau->production_vertex()->particles_end(HepMC::parents); mother++ ) {
292  if((*mother)->pdg_id() == tau->pdg_id()) return GetMother(*mother);
293  return (*mother);
294  }
295  }
296  return tau;
297 }
298 
299 
300 const std::vector<HepMC::GenParticle*> TauValidation::GetMothers(const HepMC::GenParticle* boson){
301  std::vector<HepMC::GenParticle*> mothers;
302  if ( boson->production_vertex() ) {
303  HepMC::GenVertex::particle_iterator mother;
304  for (mother = boson->production_vertex()->particles_begin(HepMC::parents); mother!= boson->production_vertex()->particles_end(HepMC::parents); mother++ ) {
305  if((*mother)->pdg_id() == boson->pdg_id()) return GetMothers(*mother);
306  mothers.push_back(*mother);
307  }
308  }
309  return mothers;
310 }
311 
312 
314  int mother_pid = 0;
315  if ( tau->production_vertex() ) {
316  HepMC::GenVertex::particle_iterator mother;
317  for (mother = tau->production_vertex()->particles_begin(HepMC::parents); mother!= tau->production_vertex()->particles_end(HepMC::parents); mother++ ) {
318  mother_pid = (*mother)->pdg_id();
319  if(mother_pid == tau->pdg_id()) return findMother(*mother); //mother_pid = -1; Make recursive to look for last tau in chain
320  }
321  }
322  return mother_pid;
323 }
324 
325 
327  if ( tau->end_vertex() ) {
328  HepMC::GenVertex::particle_iterator dau;
329  for (dau = tau->end_vertex()->particles_begin(HepMC::children); dau!= tau->end_vertex()->particles_end(HepMC::children); dau++ ) {
330  int dau_pid = (*dau)->pdg_id();
331  if(dau_pid == tau->pdg_id()) return false;
332  }
333  }
334  return true;
335 }
336 
337 
338 void TauValidation::findTauList(const HepMC::GenParticle* tau,std::vector<const HepMC::GenParticle*> &TauList){
339  TauList.insert(TauList.begin(),tau);
340  if ( tau->production_vertex() ) {
341  HepMC::GenVertex::particle_iterator mother;
342  for (mother = tau->production_vertex()->particles_begin(HepMC::parents); mother!= tau->production_vertex()->particles_end(HepMC::parents);mother++) {
343  if((*mother)->pdg_id() == tau->pdg_id()){
344  findTauList(*mother,TauList);
345  }
346  }
347  }
348 }
349 
350 void TauValidation::findFSRandBrem(const HepMC::GenParticle* p, bool doBrem, std::vector<const HepMC::GenParticle*> &ListofFSR,
351  std::vector<const HepMC::GenParticle*> &ListofBrem){
352  // note this code split the FSR and Brem based one if the tau decays into a tau+photon or not with the Fortran Tauola Interface, this is not 100% correct because photos puts the tau with the regular tau decay products.
353  if(abs(p->pdg_id())==15){
354  if(isLastTauinChain(p)){ doBrem=true;}
355  else{ doBrem=false;}
356  }
357  int photo_ID=22;
358  if ( p->end_vertex() ) {
359  HepMC::GenVertex::particle_iterator dau;
360  for (dau = p->end_vertex()->particles_begin(HepMC::children); dau!= p->end_vertex()->particles_end(HepMC::children); dau++ ) {
361  //if(doBrem) std::cout << "true " << (*dau)->pdg_id() << " " << std::endl;
362  //else std::cout << "false " << (*dau)->pdg_id() << " " << std::endl;
363  if(abs((*dau)->pdg_id()) == abs(photo_ID) && !doBrem){ListofFSR.push_back(*dau);}
364  if(abs((*dau)->pdg_id()) == abs(photo_ID) && doBrem){ListofBrem.push_back(*dau);}
365  if((*dau)->end_vertex() && (*dau)->end_vertex()->particles_out_size()>0 && abs((*dau)->pdg_id()) != 111 && abs((*dau)->pdg_id()) != 221/* remove pi0 and eta decays*/){
366  findFSRandBrem(*dau,doBrem,ListofFSR,ListofBrem);
367  }
368  }
369  }
370 }
371 
372 
373 
374 void TauValidation::FindPhotosFSR(const HepMC::GenParticle* p,std::vector<const HepMC::GenParticle*> &ListofFSR,double &BosonScale){
375  BosonScale=0.0;
376  const HepMC::GenParticle* m=GetMother(p);
377  double mother_pid=m->pdg_id();
378  if(m->end_vertex() && mother_pid!=p->pdg_id()){
379  HepMC::GenVertex::particle_iterator dau;
380  for (dau = m->end_vertex()->particles_begin(HepMC::children); dau!= m->end_vertex()->particles_end(HepMC::children); dau++ ) {
381  int dau_pid = (*dau)->pdg_id();
382  if(fabs(dau_pid) == 22) {
383  ListofFSR.push_back(*dau);
384  }
385  }
386  }
387  if(abs(mother_pid) == 24) BosonScale=1.0; // W
388  if(abs(mother_pid) == 23) BosonScale=2.0; // Z;
389  if(abs(mother_pid) == 22) BosonScale=2.0; // gamma;
390  if(abs(mother_pid) == 25) BosonScale=2.0; // HSM;
391  if(abs(mother_pid) == 35) BosonScale=2.0; // H0;
392  if(abs(mother_pid) == 36) BosonScale=2.0; // A0;
393  if(abs(mother_pid) == 37) BosonScale=1.0; //Hpm;
394 }
395 
396 
398 
399  if(abs(tau->pdg_id()) != 15 ) return -3;
400 
401  int mother_pid = findMother(tau);
402  if(mother_pid == -2) return -2;
403 
404  int label = other;
405  if(abs(mother_pid) == 24) label = W;
406  if(abs(mother_pid) == 23) label = Z;
407  if(abs(mother_pid) == 22) label = gamma;
408  if(abs(mother_pid) == 25) label = HSM;
409  if(abs(mother_pid) == 35) label = H0;
410  if(abs(mother_pid) == 36) label = A0;
411  if(abs(mother_pid) == 37) label = Hpm;
412 
413  int mother_shortpid=(abs(mother_pid)%10000);
414  if(mother_shortpid>500 && mother_shortpid<600 )label = B;
415  if(mother_shortpid>400 && mother_shortpid<500)label = D;
416  TauMothers->Fill(label,weight);
417  if(label==B || label == D || label == other) return -1;
418 
419  return mother_pid;
420 }
421 
423  int nProngs = 0;
424  if ( tau->end_vertex() ) {
425  HepMC::GenVertex::particle_iterator des;
426  for(des = tau->end_vertex()->particles_begin(HepMC::descendants);
427  des!= tau->end_vertex()->particles_end(HepMC::descendants);++des ) {
428  int pid = (*des)->pdg_id();
429  if(abs(pid) == 15) return tauProngs(*des, weight);
430  if((*des)->status() != 1) continue; // dont count unstable particles
431 
432  const HepPDT::ParticleData* pd = fPDGTable->particle((*des)->pdg_id ());
433  int charge = (int) pd->charge();
434  if(charge == 0) continue;
435  nProngs++;
436  }
437  }
438  TauProngs->Fill(nProngs,weight);
439  return nProngs;
440 }
441 
443 
444  int channel = undetermined;
445  if(tau->status() == 1) channel = stable;
446  int allCount = 0,
447  eCount = 0,
448  muCount = 0,
449  pi0Count = 0,
450  piCount = 0,
451  rhoCount = 0,
452  a1Count = 0,
453  KCount = 0,
454  KstarCount = 0;
455 
456  if ( tau->end_vertex() ) {
457  HepMC::GenVertex::particle_iterator des;
458  for(des = tau->end_vertex()->particles_begin(HepMC::descendants);
459  des!= tau->end_vertex()->particles_end(HepMC::descendants);++des ) {
460  int pid = (*des)->pdg_id();
461  if(abs(pid) == 15) return tauDecayChannel(*des,weight);
462 
463  allCount++;
464  if(abs(pid) == 11) eCount++;
465  if(abs(pid) == 13) muCount++;
466  if(abs(pid) == 111) pi0Count++;
467  if(abs(pid) == 211) piCount++;
468  if(abs(pid) == 213) rhoCount++;
469  if(abs(pid) == 20213) a1Count++;
470  if(abs(pid) == 321) KCount++;
471  if(abs(pid) == 323) KstarCount++;
472 
473  }
474  }
475  // resonances
476  if(KCount >= 1) channel = K;
477  if(KstarCount >= 1) channel = Kstar;
478  if(a1Count >= 1) channel = a1;
479  if(rhoCount >= 1) channel = rho;
480  if(channel!=undetermined && weight!=0.0) TauDecayChannels->Fill(channel,weight);
481 
482  // final state products
483  if(piCount == 1 && pi0Count == 0) channel = pi;
484  if(piCount == 1 && pi0Count == 1) channel = pi1pi0;
485  if(piCount == 1 && pi0Count > 1) channel = pinpi0;
486  if(piCount == 3 && pi0Count == 0) channel = tripi;
487  if(piCount == 3 && pi0Count > 0) channel = tripinpi0;
488  if(eCount == 1) channel = electron;
489  if(muCount == 1) channel = muon;
490  if(weight!=0.0) TauDecayChannels->Fill(channel,weight);
491  return channel;
492 }
493 
494 
495 void TauValidation::rtau(const HepMC::GenParticle* tau,int mother, int decay, double weight){
496 
497  if(decay != pi1pi0) return; // polarization only for 1-prong hadronic taus with one neutral pion to make a clean case
498 
499  if(tau->momentum().perp() < tauEtCut) return; // rtau visible only for boosted taus
500 
501  double rTau = 0;
502  double ltrack = leadingPionMomentum(tau, weight);
503  double visibleTauE = visibleTauEnergy(tau);
504 
505  if(visibleTauE != 0) rTau = ltrack/visibleTauE;
506 
507  if(abs(mother) == 24) TauRtauW->Fill(rTau,weight);
508  if(abs(mother) == 37) TauRtauHpm->Fill(rTau,weight);
509 }
510 
511  void TauValidation::spinEffectsWHpm(const HepMC::GenParticle* tau,int mother, int decay, std::vector<HepMC::GenParticle*> &part,double weight){
512  if(decay == TauDecay::JAK_PION || decay == TauDecay::JAK_MUON || decay == TauDecay::JAK_ELECTRON){ // polarization only for 1-prong hadronic taus with no neutral pions
513  TLorentzVector momP4 = motherP4(tau);
514  TLorentzVector pionP4 = leadingPionP4(tau);
515  pionP4.Boost(-1*momP4.BoostVector());
516  double energy = pionP4.E()/(momP4.M()/2);
517  if(decay == TauDecay::JAK_PION){
518  if(abs(mother) == 24) TauSpinEffectsW_X->Fill(energy,weight);
519  if(abs(mother) == 37) TauSpinEffectsHpm_X->Fill(energy,weight);
520  }
521  if(decay == TauDecay::JAK_MUON){
522  if(abs(mother) == 24) TauSpinEffectsW_muX->Fill(energy,weight);
523  if(abs(mother) == 37) TauSpinEffectsHpm_muX->Fill(energy,weight);
524  }
525  if(decay == TauDecay::JAK_ELECTRON){
526  if(abs(mother) == 24) TauSpinEffectsW_eX->Fill(energy,weight);
527  if(abs(mother) == 37) TauSpinEffectsHpm_eX->Fill(energy,weight);
528  }
529 
530  }
531  else if(decay==TauDecay::JAK_RHO_PIPI0){
532  TLorentzVector rho(0,0,0,0),pi(0,0,0,0);
533  for(unsigned int i=0;i<part.size();i++){
534  TLorentzVector LV(part.at(i)->momentum().px(),part.at(i)->momentum().py(),part.at(i)->momentum().pz(),part.at(i)->momentum().e());
535  if(abs(part.at(i)->pdg_id())==PdtPdgMini::pi_plus){pi+=LV; rho+=LV;}
536  if(abs(part.at(i)->pdg_id())==PdtPdgMini::pi0){rho+=LV;}
537  }
538  if(abs(mother) == 24) TauSpinEffectsW_UpsilonRho->Fill(2*pi.Pt()/rho.Pt()-1,weight);
539  if(abs(mother) == 37) TauSpinEffectsHpm_UpsilonRho->Fill(2*pi.Pt()/rho.Pt()-1,weight);
540  }
541  else if(decay==TauDecay::JAK_A1_3PI){ // only for pi2pi0 for now
542  TLorentzVector a1(0,0,0,0),pi_p(0,0,0,0),pi_m(0,0,0,0);
543  int nplus(0),nminus(0);
544  for(unsigned int i=0;i<part.size();i++){
545  TLorentzVector LV(part.at(i)->momentum().px(),part.at(i)->momentum().py(),part.at(i)->momentum().pz(),part.at(i)->momentum().e());
546  if(part.at(i)->pdg_id()==PdtPdgMini::pi_plus){ pi_p+=LV; a1+=LV; nplus++;}
547  if(part.at(i)->pdg_id()==PdtPdgMini::pi_minus){pi_m+=LV; a1+=LV; nminus++;}
548  }
549  double gamma=0;
550  if(nplus+nminus==3 && nplus==1) gamma=2*pi_p.Pt()/a1.Pt()-1;
551  if(nplus+nminus==3 && nminus==1) gamma=2*pi_m.Pt()/a1.Pt()-1;
552  else{
553  pi_p+=pi_m; gamma=2*pi_p.Pt()/a1.Pt()-1;
554  }
555  if(abs(mother) == 24) TauSpinEffectsW_UpsilonA1->Fill(gamma,weight);
556  if(abs(mother) == 37) TauSpinEffectsHpm_UpsilonA1->Fill(gamma,weight);
557  }
558 }
559 
561  TLorentzVector tautau(0,0,0,0);
562  TLorentzVector pipi(0,0,0,0);
563  TLorentzVector taum(0,0,0,0);
564  TLorentzVector taup(0,0,0,0);
565  TLorentzVector rho_plus,rho_minus,pi_rhominus,pi0_rhominus,pi_rhoplus,pi0_rhoplus,pi_plus,pi_minus;
566  bool hasrho_minus(false),hasrho_plus(false),haspi_minus(false),haspi_plus(false);
567  int nSinglePionDecays(0),nSingleMuonDecays(0),nSingleElectronDecays(0);
568  double x1(0),x2(0);
569  TLorentzVector Zboson(boson->momentum().px(),boson->momentum().py(),boson->momentum().pz(),boson->momentum().e());
570  if ( boson->end_vertex() ) {
571  HepMC::GenVertex::particle_iterator des;
572  for(des = boson->end_vertex()->particles_begin(HepMC::children); des!= boson->end_vertex()->particles_end(HepMC::children);++des ) {
573  int pid = (*des)->pdg_id();
574  if(abs(findMother(*des)) != 15 && abs(pid) == 15 && (tauDecayChannel(*des) == pi || tauDecayChannel(*des) == muon || tauDecayChannel(*des) == electron )){
575  if(tauDecayChannel(*des) == pi)nSinglePionDecays++;
576  if(tauDecayChannel(*des) == muon)nSingleMuonDecays++;
577  if(tauDecayChannel(*des) == electron)nSingleElectronDecays++;
578  TLorentzVector LVtau((*des)->momentum().px(),(*des)->momentum().py(),(*des)->momentum().pz(),(*des)->momentum().e());
579  tautau += LVtau;
580  TLorentzVector LVpi=leadingPionP4(*des);
581  pipi+=LVpi;
582  const HepPDT::ParticleData* pd = fPDGTable->particle((*des)->pdg_id ());
583  int charge = (int) pd->charge();
584  LVtau.Boost(-1*Zboson.BoostVector());
585  LVpi.Boost(-1*Zboson.BoostVector());
586  if(tauDecayChannel(*des) == pi){
587  if(abs(boson->pdg_id())==PdtPdgMini::Z0) TauSpinEffectsZ_X->Fill(LVpi.P()/LVtau.E(),weight);
588  if(abs(boson->pdg_id())==PdtPdgMini::Higgs0) TauSpinEffectsH_X->Fill(LVpi.P()/LVtau.E(),weight);
589  }
590  if(charge<0){x1=LVpi.P()/LVtau.E(); taum=LVtau;}
591  else{ x2=LVpi.P()/LVtau.E();}
592  }
593  if(abs(findMother(*des)) != 15 && abs(pid) == 15 && (tauDecayChannel(*des) == rho || tauDecayChannel(*des) == pi1pi0) ){
594  if ( (*des)->end_vertex() ) {
595  HepMC::GenVertex::particle_iterator tauprod;
596  TLorentzVector LVtau((*des)->momentum().px(),(*des)->momentum().py(),(*des)->momentum().pz(),(*des)->momentum().e());
597  if(pid == 15)taum=LVtau;
598  if(pid ==-15)taup=LVtau;
599  for(tauprod = (*des)->end_vertex()->particles_begin(HepMC::descendants); tauprod!= (*des)->end_vertex()->particles_end(HepMC::descendants);++tauprod ) {
600  int pid_d = (*tauprod)->pdg_id();
601  if(abs(pid_d)==211 || abs(pid_d)==111){
602  TLorentzVector LV((*tauprod)->momentum().px(),(*tauprod)->momentum().py(),(*tauprod)->momentum().pz(),(*tauprod)->momentum().e());
603  if(pid==15){
604  hasrho_minus=true;
605  if(pid_d==-211 ){ pi_rhominus=LV;}
606  if(abs(pid_d)==111 ){ pi0_rhominus=LV;}
607  }
608  if(pid==-15){
609  hasrho_plus=true;
610  if(pid_d==211 ){pi_rhoplus=LV;}
611  if(abs(pid_d)==111 ){pi0_rhoplus=LV;}
612  }
613  }
614  }
615  }
616  }
617  if(abs(findMother(*des)) != 15 && abs(pid) == 15 && (tauDecayChannel(*des) == pi) ){
618  if ( (*des)->end_vertex() ) {
619  HepMC::GenVertex::particle_iterator tauprod;
620  TLorentzVector LVtau((*des)->momentum().px(),(*des)->momentum().py(),(*des)->momentum().pz(),(*des)->momentum().e());
621  if(pid == 15)taum=LVtau;
622  if(pid ==-15)taup=LVtau;
623  for(tauprod = (*des)->end_vertex()->particles_begin(HepMC::descendants); tauprod!= (*des)->end_vertex()->particles_end(HepMC::descendants);++tauprod ) {
624  int pid_d = (*tauprod)->pdg_id();
625  if(abs(pid_d)==211 || abs(pid_d)==111){
626  TLorentzVector LV((*tauprod)->momentum().px(),(*tauprod)->momentum().py(),(*tauprod)->momentum().pz(),(*tauprod)->momentum().e());
627  if(pid==15){
628  haspi_minus=true;
629  if(pid_d==-211 ){ pi_minus=LV;}
630  }
631  if(pid==-15){
632  haspi_plus=true;
633  if(pid_d==211 ){pi_plus=LV;}
634  }
635  }
636  }
637  }
638  }
639  }
640  }
641  if(hasrho_minus && hasrho_plus){
642  //compute rhorho
643  rho_minus=pi_rhominus;
644  rho_minus+=pi0_rhominus;
645  rho_plus=pi_rhoplus;
646  rho_plus+=pi0_rhoplus;
647  TLorentzVector rhorho=rho_minus;rhorho+=rho_plus;
648 
649  // boost to rhorho cm
650  TLorentzVector pi_rhoplusb=pi_rhoplus; pi_rhoplusb.Boost(-1*rhorho.BoostVector());
651  TLorentzVector pi0_rhoplusb=pi0_rhoplus; pi0_rhoplusb.Boost(-1*rhorho.BoostVector());
652  TLorentzVector pi_rhominusb=pi_rhominus; pi_rhominusb.Boost(-1*rhorho.BoostVector());
653  TLorentzVector pi0_rhominusb=pi0_rhominus; pi0_rhominusb.Boost(-1*rhorho.BoostVector());
654 
655  // compute n+/-
656  TVector3 n_plus=pi_rhoplusb.Vect().Cross(pi0_rhoplusb.Vect());
657  TVector3 n_minus=pi_rhominusb.Vect().Cross(pi0_rhominusb.Vect());
658 
659  // compute the acoplanarity
660  double Acoplanarity=acos(n_plus.Dot(n_minus)/(n_plus.Mag()*n_minus.Mag()));
661  if(pi_rhominus.Vect().Dot(n_plus)>0){Acoplanarity*=-1;Acoplanarity+=2*TMath::Pi();}
662 
663  // now boost to tau frame
664  pi_rhoplus.Boost(-1*taup.BoostVector());
665  pi0_rhoplus.Boost(-1*taup.BoostVector());
666  pi_rhominus.Boost(-1*taum.BoostVector());
667  pi0_rhominus.Boost(-1*taum.BoostVector());
668 
669  // compute y1 and y2
670  double y1=(pi_rhoplus.E()-pi0_rhoplus.E())/(pi_rhoplus.E()+pi0_rhoplus.E());
671  double y2=(pi_rhominus.E()-pi0_rhominus.E())/(pi_rhominus.E()+pi0_rhominus.E());
672 
673  // fill histograms
674  if(abs(boson->pdg_id())==PdtPdgMini::Higgs0 && y1*y2<0) TauSpinEffectsH_rhorhoAcoplanarityminus->Fill(Acoplanarity,weight);
675  if(abs(boson->pdg_id())==PdtPdgMini::Higgs0 && y1*y2>0) TauSpinEffectsH_rhorhoAcoplanarityplus->Fill(Acoplanarity,weight);
676  }
677  if(haspi_minus && haspi_plus){
678  TLorentzVector tauporig=taup;
679  TLorentzVector taumorig=taum;
680 
681  // now boost to Higgs frame
682  pi_plus.Boost(-1*Zboson.BoostVector());
683  pi_minus.Boost(-1*Zboson.BoostVector());
684 
685  taup.Boost(-1*Zboson.BoostVector());
686  taum.Boost(-1*Zboson.BoostVector());
687 
688  if(abs(boson->pdg_id())==PdtPdgMini::Higgs0){
689  TauSpinEffectsH_pipiAcollinearity->Fill(acos(pi_plus.Vect().Dot(pi_minus.Vect())/(pi_plus.P()*pi_minus.P())));
690  TauSpinEffectsH_pipiAcollinearityzoom->Fill(acos(pi_plus.Vect().Dot(pi_minus.Vect())/(pi_plus.P()*pi_minus.P())));
691  }
692 
693  double proj_m=taum.Vect().Dot(pi_minus.Vect())/(taum.P()*taup.P());
694  double proj_p=taup.Vect().Dot(pi_plus.Vect())/(taup.P()*taup.P());
695  TVector3 Tau_m=taum.Vect();
696  TVector3 Tau_p=taup.Vect();
697  Tau_m*=proj_m;
698  Tau_p*=proj_p;
699  TVector3 Pit_m=pi_minus.Vect()-Tau_m;
700  TVector3 Pit_p=pi_plus.Vect()-Tau_p;
701 
702  double Acoplanarity=acos(Pit_m.Dot(Pit_p)/(Pit_p.Mag()*Pit_m.Mag()));
703  TVector3 n=Pit_p.Cross(Pit_m);
704  if(n.Dot(Tau_m)/Tau_m.Mag()>0){Acoplanarity*=-1; Acoplanarity+=2*TMath::Pi();}
705  // fill histograms
706  if(abs(boson->pdg_id())==PdtPdgMini::Higgs0) TauSpinEffectsH_pipiAcoplanarity->Fill(Acoplanarity,weight);
707  taup=tauporig;
708  taum=taumorig;
709 
710  }
711 
712 
713  if(nSingleMuonDecays==2){
714  if(abs(boson->pdg_id())==PdtPdgMini::Z0) TauSpinEffectsZ_muX->Fill(x1,weight);
715  if(abs(boson->pdg_id())==PdtPdgMini::Higgs0) TauSpinEffectsH_muX->Fill(x1,weight);
716  }
717  if(nSingleElectronDecays==2){
718  if(abs(boson->pdg_id())==PdtPdgMini::Z0) TauSpinEffectsZ_eX->Fill(x1,weight);
719  if(abs(boson->pdg_id())==PdtPdgMini::Higgs0) TauSpinEffectsH_eX->Fill(x1,weight);
720  }
721  if(nSinglePionDecays == 2 && tautau.M()!= 0) {
722  for(int i=0;i<zsbins;i++){
723  double zslow=((double)i)*(zsmax-zsmin)/((double)zsbins)+zsmin;
724  double zsup=((double)i+1)*(zsmax-zsmin)/((double)zsbins)+zsmin;
725  double aup=Zstoa(zsup), alow=Zstoa(zslow);
726  if(x2-x1>alow && x2-x1<aup){
727  double zs=(zsup+zslow)/2;
728  if(abs(boson->pdg_id())==PdtPdgMini::Z0) TauSpinEffectsZ_Zs->Fill(zs,weight);
729  if(abs(boson->pdg_id())==PdtPdgMini::Higgs0) TauSpinEffectsH_Zs->Fill(zs,weight);
730  break;
731  }
732  }
733  if(abs(boson->pdg_id())==PdtPdgMini::Z0) TauSpinEffectsZ_MVis->Fill(pipi.M()/tautau.M(),weight);
734  if(abs(boson->pdg_id())==PdtPdgMini::Higgs0) TauSpinEffectsH_MVis->Fill(pipi.M()/tautau.M(),weight);
735 
736  if(x1!=0){
737  const std::vector<HepMC::GenParticle*> m=GetMothers(boson);
738  int q(0),qbar(0);
739  TLorentzVector Z(0,0,0,0);
740  for(unsigned int i=0;i<m.size();i++){
741  if(m.at(i)->pdg_id()==PdtPdgMini::d || m.at(i)->pdg_id()==PdtPdgMini::u ){q++;}
742  if(m.at(i)->pdg_id()==PdtPdgMini::anti_d || m.at(i)->pdg_id()==PdtPdgMini::anti_u ){qbar++;}
743  }
744  if(q==1 && qbar==1){// assume q has largest E (valence vs see quarks)
745  if(taum.Vect().Dot(Zboson.Vect())/(Zboson.P()*taum.P())>0){
746  if(abs(boson->pdg_id())==PdtPdgMini::Z0) TauSpinEffectsZ_Xf->Fill(x1,weight);
747  if(abs(boson->pdg_id())==PdtPdgMini::Higgs0) TauSpinEffectsH_Xf->Fill(x1,weight);
748  }
749  else{
750  if(abs(boson->pdg_id())==PdtPdgMini::Z0) TauSpinEffectsZ_Xb->Fill(x1,weight);
751  if(abs(boson->pdg_id())==PdtPdgMini::Higgs0) TauSpinEffectsH_Xb->Fill(x1,weight);
752  }
753  }
754  }
755  }
756 }
757 
758 double TauValidation::Zstoa(double zs){
759  double a=1-sqrt(fabs(1.0-2*fabs(zs)));
760  if(zs<0){
761  a*=-1.0;
762  }
763  return a;
764 }
765 
766 
768  return leadingPionP4(tau).P();
769 }
770 
772 
773  TLorentzVector p4(0,0,0,0);
774 
775  if ( tau->end_vertex() ) {
776  HepMC::GenVertex::particle_iterator des;
777  for(des = tau->end_vertex()->particles_begin(HepMC::descendants);
778  des!= tau->end_vertex()->particles_end(HepMC::descendants);++des ) {
779  int pid = (*des)->pdg_id();
780 
781  if(abs(pid) == 15) return leadingPionP4(*des);
782 
783  if(!(abs(pid)==211 || abs(pid)==13 || abs(pid)==11)) continue;
784 
785  if((*des)->momentum().rho() > p4.P()) {
786  p4 = TLorentzVector((*des)->momentum().px(),
787  (*des)->momentum().py(),
788  (*des)->momentum().pz(),
789  (*des)->momentum().e());
790  }
791  }
792  }
793 
794  return p4;
795 }
796 
798  const HepMC::GenParticle* m=GetMother(tau);
799  return TLorentzVector(m->momentum().px(),m->momentum().py(),m->momentum().pz(),m->momentum().e());
800 }
801 
803  TLorentzVector p4(tau->momentum().px(),
804  tau->momentum().py(),
805  tau->momentum().pz(),
806  tau->momentum().e());
807 
808  if ( tau->end_vertex() ) {
809  HepMC::GenVertex::particle_iterator des;
810  for(des = tau->end_vertex()->particles_begin(HepMC::descendants);
811  des!= tau->end_vertex()->particles_end(HepMC::descendants);++des ) {
812  int pid = (*des)->pdg_id();
813 
814  if(abs(pid) == 15) return visibleTauEnergy(*des);
815 
816  if(abs(pid) == 12 || abs(pid) == 14 || abs(pid) == 16) {
817  p4 -= TLorentzVector((*des)->momentum().px(),
818  (*des)->momentum().py(),
819  (*des)->momentum().pz(),
820  (*des)->momentum().e());
821  }
822  }
823  }
824 
825  return p4.E();
826 }
827 
829  // Find First tau in chain
830  std::vector<const HepMC::GenParticle*> TauList;
831  findTauList(tau,TauList);
832 
833  // Get List of Gauge Boson to tau(s) FSR and Brem
834  bool passedW=false;
835  std::vector<const HepMC::GenParticle*> ListofFSR; ListofFSR.clear();
836  std::vector<const HepMC::GenParticle*> ListofBrem; ListofBrem.clear();
837  std::vector<const HepMC::GenParticle*> FSR_photos; FSR_photos.clear();
838  double BosonScale(1);
839  if(TauList.size()>0){
840  TauValidation::findFSRandBrem(TauList.at(0),passedW,ListofFSR,ListofBrem);
841  TauValidation::FindPhotosFSR(TauList.at(0),FSR_photos,BosonScale);
842 
843  // Add the Tau Brem. information
844  TauBremPhotonsN->Fill(ListofBrem.size(),weight);
845  double photonPtSum=0;
846  for(unsigned int i=0;i<ListofBrem.size();i++){
847  photonPtSum+=ListofBrem.at(i)->momentum().perp();
848  TauBremPhotonsPt->Fill(ListofBrem.at(i)->momentum().perp(),weight);
849  }
850  TauBremPhotonsPtSum->Fill(photonPtSum,weight);
851 
852  // Now add the Gauge Boson FSR information
853  if(BosonScale!=0){
854  TauFSRPhotonsN->Fill(ListofFSR.size(),weight);
855  photonPtSum=0;
856  for(unsigned int i=0;i<ListofFSR.size();i++){
857  photonPtSum+=ListofFSR.at(i)->momentum().perp();
858  TauFSRPhotonsPt->Fill(ListofFSR.at(i)->momentum().perp(),weight);
859  }
860  double FSR_photosSum(0);
861  for(unsigned int i=0;i<FSR_photos.size();i++){
862  FSR_photosSum+=FSR_photos.at(i)->momentum().perp();
863  TauFSRPhotonsPt->Fill(FSR_photos.at(i)->momentum().perp()/BosonScale,weight*BosonScale);
864  }
865  TauFSRPhotonsPtSum->Fill(photonPtSum+FSR_photosSum/BosonScale,weight);
866  }
867  }
868 }
869 
void findTauList(const HepMC::GenParticle *tau, std::vector< const HepMC::GenParticle * > &TauList)
MonitorElement * TauFSRPhotonsN
const double beta
const double Pi
int i
Definition: DBlmapReader.cc:9
virtual ~TauValidation()
MonitorElement * TauSpinEffectsH_rhorhoAcoplanarityplus
MonitorElement * TauPhi
MonitorElement * TauSpinEffectsW_X
TPRegexp parents
Definition: eve_filter.cc:24
MonitorElement * TauSpinEffectsHpm_UpsilonA1
std::vector< std::vector< MonitorElement * > > JAKInvMass
MonitorElement * TauSpinEffectsZ_Xb
MonitorElement * nTaus
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:872
WeightManager wmanager_
Definition: TauValidation.h:73
MonitorElement * DecayLength
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
MonitorElement * TauRtauW
TauValidation(const edm::ParameterSet &)
MonitorElement * TauPt
MonitorElement * TauSpinEffectsH_MVis
void spinEffectsZH(const HepMC::GenParticle *boson, double weight)
MonitorElement * TauMothers
MonitorElement * TauSpinEffectsHpm_eX
MonitorElement * TauSpinEffectsW_muX
MonitorElement * TauBremPhotonsN
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
MonitorElement * JAKID
MonitorElement * TauSpinEffectsZ_Xf
int tauProngs(const HepMC::GenParticle *, double weight)
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
MonitorElement * TauBremPhotonsPtSum
virtual void endJob()
MonitorElement * TauSpinEffectsH_pipiAcoplanarity
MonitorElement * TauSpinEffectsHpm_UpsilonRho
double charge(const std::vector< uint8_t > &Ampls)
virtual void analyze(const edm::Event &, const edm::EventSetup &)
void getData(T &iHolder) const
Definition: EventSetup.h:67
MonitorElement * TauSpinEffectsZ_eX
int findMother(const HepMC::GenParticle *)
void Fill(long long x)
MonitorElement * TauProngs
int tauMother(const HepMC::GenParticle *, double weight)
int tauDecayChannel(const HepMC::GenParticle *, double weight=0.0)
bool isTauFinalStateParticle(int pdgid)
Definition: TauDecay.cc:38
MonitorElement * TauSpinEffectsH_pipiAcollinearityzoom
const std::vector< HepMC::GenParticle * > GetMothers(const HepMC::GenParticle *boson)
math::XYZTLorentzVectorD LV
int iEvent
Definition: GenABIO.cc:243
MonitorElement * TauSpinEffectsW_UpsilonRho
void photons(const HepMC::GenParticle *, double weight)
MonitorElement * TauSpinEffectsH_rhorhoAcoplanarityminus
T sqrt(T t)
Definition: SSEVec.h:48
double p4[4]
Definition: TauolaWrapper.h:92
double leadingPionMomentum(const HepMC::GenParticle *, double weight)
virtual void beginJob()
MonitorElement * TauSpinEffectsH_Xf
MonitorElement * TauFSRPhotonsPt
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
MonitorElement * TauSpinEffectsZ_muX
MonitorElement * TauDecayChannels
edm::InputTag hepmcCollection_
Definition: TauValidation.h:96
MonitorElement * TauSpinEffectsH_Zs
MonitorElement * TauFSRPhotonsPtSum
HepPDT::ParticleData ParticleData
double visibleTauEnergy(const HepMC::GenParticle *)
TLorentzVector leadingPionP4(const HepMC::GenParticle *)
bool isLastTauinChain(const HepMC::GenParticle *tau)
MonitorElement * TauSpinEffectsZ_MVis
MonitorElement * TauSpinEffectsH_pipiAcollinearity
TLorentzVector motherP4(const HepMC::GenParticle *)
void rtau(const HepMC::GenParticle *, int, int, double weight)
part
Definition: HCALResponse.h:20
tuple pid
Definition: sysUtil.py:22
void FindPhotosFSR(const HepMC::GenParticle *p, std::vector< const HepMC::GenParticle * > &ListofFSR, double &BosonScale)
MonitorElement * TauSpinEffectsZ_X
virtual void endRun(const edm::Run &, const edm::EventSetup &)
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.
MonitorElement * TauSpinEffectsH_eX
MonitorElement * TauSpinEffectsW_eX
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
double Zstoa(double zs)
const HepMC::GenParticle * GetMother(const HepMC::GenParticle *tau)
void spinEffectsWHpm(const HepMC::GenParticle *, int, int, std::vector< HepMC::GenParticle * > &part, double weight)
MonitorElement * LifeTime
MonitorElement * TauSpinEffectsH_muX
MonitorElement * TauSpinEffectsW_UpsilonA1
double a
Definition: hdecay.h:121
MonitorElement * TauEta
DQMStore * dbe
ME&#39;s &quot;container&quot;.
void findFSRandBrem(const HepMC::GenParticle *p, bool doBrem, std::vector< const HepMC::GenParticle * > &ListofFSR, std::vector< const HepMC::GenParticle * > &ListofBrem)
unsigned int NJAKID
MonitorElement * TauSpinEffectsZ_Zs
MonitorElement * TauSpinEffectsHpm_X
virtual void beginRun(const edm::Run &, const edm::EventSetup &)
int weight
Definition: histoStyle.py:50
MonitorElement * TauSpinEffectsHpm_muX
double weight(const edm::Event &)
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
MonitorElement * TauSpinEffectsH_Xb
std::vector< HepMC::GenParticle * > Get_TauDecayProducts()
bool AnalyzeTau(HepMC::GenParticle *Tau, unsigned int &JAK_ID, unsigned int &TauBitMask, bool dores=true, bool dopi0=true)
MonitorElement * nPrimeTaus
MonitorElement * TauBremPhotonsPt
MonitorElement * TauRtauHpm
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:584
Definition: Run.h:41
unsigned int nProng(unsigned int &TauBitMask)
Definition: TauDecay.h:103
MonitorElement * TauSpinEffectsH_X