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