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 #include "CLHEP/Units/defs.h"
8 #include "CLHEP/Units/PhysicalConstants.h"
15 #include <iostream>
16 using namespace edm;
17 
19  wmanager_(iPSet,consumesCollector())
20  ,genparticleCollection_(iPSet.getParameter<edm::InputTag>("genparticleCollection"))
21  ,hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection"))
22  ,NJAKID(22)
23  ,zsbins(20)
24  ,zsmin(-0.5)
25  ,zsmax(0.5)
26 {
27  hepmcCollectionToken_=consumes<HepMCProduct>(hepmcCollection_);
28  genparticleCollectionToken_=consumes<reco::GenParticleCollection>(genparticleCollection_);
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);
66 
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 = i.book1D("DecayLength","#tau Decay Length", 100 ,-20,20); DecayLength->setAxisTitle("L_{#tau} (cm)");
79  LifeTime = i.book1D("LifeTime","#tau LifeTime ", 500 ,0,10000E-15); LifeTime->setAxisTitle("#tau_{#tau} (s)");
80 
81  TauSpinEffectsW_X = i.book1D("TauSpinEffectsWX","Pion energy in W rest frame", 50 ,0,1); TauSpinEffectsW_X->setAxisTitle("X");
82  TauSpinEffectsHpm_X = i.book1D("TauSpinEffectsHpmX","Pion energy in Hpm rest frame", 50 ,0,1); TauSpinEffectsHpm_X->setAxisTitle("X");
83 
84  TauSpinEffectsW_eX = i.book1D("TauSpinEffectsWeX","e energy in W rest frame", 50 ,0,1); TauSpinEffectsW_eX->setAxisTitle("X");
85  TauSpinEffectsHpm_eX = i.book1D("TauSpinEffectsHpmeX","e energy in Hpm rest frame", 50 ,0,1); TauSpinEffectsHpm_eX->setAxisTitle("X");
86 
87  TauSpinEffectsW_muX = i.book1D("TauSpinEffectsWmuX","mu energy in W rest frame", 50 ,0,1); TauSpinEffectsW_muX->setAxisTitle("X");
88  TauSpinEffectsHpm_muX = i.book1D("TauSpinEffectsHpmmuX","mu energy in Hpm rest frame", 50 ,0,1); TauSpinEffectsHpm_muX->setAxisTitle("X");
89 
90  TauSpinEffectsW_UpsilonRho = i.book1D("TauSpinEffectsWUpsilonRho","#Upsilon for #rho", 50 ,-1,1); TauSpinEffectsW_UpsilonRho->setAxisTitle("#Upsilon");
91  TauSpinEffectsHpm_UpsilonRho = i.book1D("TauSpinEffectsHpmUpsilonRho","#Upsilon for #rho", 50 ,-1,1); TauSpinEffectsHpm_UpsilonRho->setAxisTitle("#Upsilon");
92 
93  TauSpinEffectsW_UpsilonA1 = i.book1D("TauSpinEffectsWUpsilonA1","#Upsilon for a1", 50 ,-1,1); TauSpinEffectsW_UpsilonA1->setAxisTitle("#Upsilon");
94  TauSpinEffectsHpm_UpsilonA1 = i.book1D("TauSpinEffectsHpmUpsilonA1","#Upsilon for a1", 50 ,-1,1); TauSpinEffectsHpm_UpsilonA1->setAxisTitle("#Upsilon");
95 
96  TauSpinEffectsH_pipiAcoplanarity = i.book1D("TauSpinEffectsH_pipiAcoplanarity","H Acoplanarity for #pi^{-}#pi^{+}", 50 ,0,2*TMath::Pi()); TauSpinEffectsH_pipiAcoplanarity->setAxisTitle("Acoplanarity");
97 
98  TauSpinEffectsH_pipiAcollinearity = i.book1D("TauSpinEffectsH_pipiAcollinearity","H Acollinearity for #pi^{-}#pi^{+}", 50 ,0,TMath::Pi()); TauSpinEffectsH_pipiAcollinearity->setAxisTitle("Acollinearity");
99  TauSpinEffectsH_pipiAcollinearityzoom = i.book1D("TauSpinEffectsH_pipiAcollinearityzoom","H Acollinearity for #pi^{-}#pi^{+}", 50 ,3,TMath::Pi()); TauSpinEffectsH_pipiAcollinearityzoom->setAxisTitle("Acollinearity");
100 
101  TauSpinEffectsZ_MVis = i.book1D("TauSpinEffectsZMVis","Mass of pi+ pi-", 25 ,0,1.1); TauSpinEffectsZ_MVis->setAxisTitle("M_{#pi^{+}#pi^{-}}");
102  TauSpinEffectsH_MVis = i.book1D("TauSpinEffectsHMVis","Mass of pi+ pi-", 25 ,0,1.1); TauSpinEffectsZ_MVis->setAxisTitle("M_{#pi^{+}#pi^{-}}");
103 
104  TauSpinEffectsZ_Zs = i.book1D("TauSpinEffectsZZs","Z_{s}", zsbins ,zsmin,zsmax); TauSpinEffectsZ_Zs->setAxisTitle("Z_{s}");
105  TauSpinEffectsH_Zs = i.book1D("TauSpinEffectsHZs","Z_{s}", zsbins ,zsmin,zsmax); TauSpinEffectsZ_Zs->setAxisTitle("Z_{s}");
106 
107  TauSpinEffectsZ_X= i.book1D("TauSpinEffectsZX","X of #tau^{-}", 25 ,0,1.0); TauSpinEffectsZ_X->setAxisTitle("X");
108  TauSpinEffectsZ_X50to75= i.book1D("TauSpinEffectsZX50to75","X of #tau^{-} (50GeV-75GeV)", 10 ,0,1.0); TauSpinEffectsZ_X->setAxisTitle("X");
109  TauSpinEffectsZ_X75to88= i.book1D("TauSpinEffectsZX75to88","X of #tau^{-} (75GeV-88GeV)", 10 ,0,1.0); TauSpinEffectsZ_X->setAxisTitle("X");
110  TauSpinEffectsZ_X88to100= i.book1D("TauSpinEffectsZX88to100","X of #tau^{-} (88GeV-100GeV)", 10 ,0,1.0); TauSpinEffectsZ_X->setAxisTitle("X");
111  TauSpinEffectsZ_X100to120= i.book1D("TauSpinEffectsZX100to120","X of #tau^{-} (100GeV-120GeV)", 10 ,0,1.0); TauSpinEffectsZ_X->setAxisTitle("X");
112  TauSpinEffectsZ_X120UP= i.book1D("TauSpinEffectsZX120UP","X of #tau^{-} (>120GeV)", 10 ,0,1.0); TauSpinEffectsZ_X->setAxisTitle("X");
113 
114 
115  TauSpinEffectsH_X= i.book1D("TauSpinEffectsH_X","X of #tau^{-}", 25 ,0,1.0); TauSpinEffectsH_X->setAxisTitle("X");
116 
117  TauSpinEffectsZ_Xf = i.book1D("TauSpinEffectsZXf","X of forward emitted #tau^{-}", 25 ,0,1.0); TauSpinEffectsZ_Xf->setAxisTitle("X_{f}");
118  TauSpinEffectsH_Xf = i.book1D("TauSpinEffectsHXf","X of forward emitted #tau^{-}", 25 ,0,1.0); TauSpinEffectsZ_Xf->setAxisTitle("X_{f}");
119 
120  TauSpinEffectsZ_Xb = i.book1D("TauSpinEffectsZXb","X of backward emitted #tau^{-}", 25 ,0,1.0); TauSpinEffectsZ_Xb->setAxisTitle("X_{b}");
121  TauSpinEffectsH_Xb = i.book1D("TauSpinEffectsHXb","X of backward emitted #tau^{-}", 25 ,0,1.0); TauSpinEffectsZ_Xb->setAxisTitle("X_{b}");
122 
123  TauSpinEffectsZ_eX = i.book1D("TauSpinEffectsZeX","e energy in Z rest frame", 50 ,0,1); TauSpinEffectsZ_eX->setAxisTitle("X");
124  TauSpinEffectsH_eX = i.book1D("TauSpinEffectsHeX","e energy in H rest frame", 50 ,0,1); TauSpinEffectsH_eX->setAxisTitle("X");
125 
126  TauSpinEffectsZ_muX = i.book1D("TauSpinEffectsZmuX","mu energy in z rest frame", 50 ,0,1); TauSpinEffectsZ_muX->setAxisTitle("X");
127  TauSpinEffectsH_muX = i.book1D("TauSpinEffectsHmuX","mu energy in H rest frame", 50 ,0,1); TauSpinEffectsH_muX->setAxisTitle("X");
128 
129  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)");
130  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)");
131 
132  TauFSRPhotonsN=i.book1D("TauFSRPhotonsN","FSR Photons radiating from/with tau (Gauge Boson)", 5 ,-0.5,4.5);
133  TauFSRPhotonsN->setAxisTitle("N FSR Photons radiating from/with tau");
134  TauFSRPhotonsPt=i.book1D("TauFSRPhotonsPt","Pt of FSR Photons radiating from/with tau (Gauge Boson)", 100 ,0,100);
135  TauFSRPhotonsPt->setAxisTitle("P_{t} of FSR Photons radiating from/with tau [per tau]");
136  TauFSRPhotonsPtSum=i.book1D("TauFSRPhotonsPtSum","Pt of FSR Photons radiating from/with tau (Gauge Boson)", 100 ,0,100);
137  TauFSRPhotonsPtSum->setAxisTitle("P_{t} of FSR Photons radiating from/with tau [per tau]");
138 
139  TauBremPhotonsN=i.book1D("TauBremPhotonsN","Brem. Photons radiating in tau decay", 5 ,-0.5,4.5);
140  TauBremPhotonsN->setAxisTitle("N FSR Photons radiating from/with tau");
141  TauBremPhotonsPt=i.book1D("TauBremPhotonsPt","Sum Brem Pt ", 100 ,0,100);
142  TauFSRPhotonsPt->setAxisTitle("P_{t} of Brem. Photons radiating in tau decay");
143  TauBremPhotonsPtSum =i.book1D("TauBremPhotonsPtSum","Sum of Brem Pt ", 100 ,0,100);
144  TauFSRPhotonsPtSum->setAxisTitle("Sum P_{t} of Brem. Photons radiating in tau decay");
145 
146  JAKID =i.book1D("JAKID","JAK ID",NJAKID+1,-0.5,NJAKID+0.5);
147  for(unsigned int j=0; j<NJAKID+1;j++){
148  JAKInvMass.push_back(std::vector<MonitorElement *>());
149  TString tmp="JAKID";
150  tmp+=j;
151  JAKInvMass.at(j).push_back(i.book1D("M"+tmp,"M_{"+tmp+"} (GeV)", 80 ,0,2.0));
152  if(j==TauDecay::JAK_A1_3PI ||
155  JAKInvMass.at(j).push_back(i.book1D("M13"+tmp,"M_{13,"+tmp+"} (GeV)", 80 ,0,2.0));
156  JAKInvMass.at(j).push_back(i.book1D("M23"+tmp,"M_{23,"+tmp+"} (GeV)", 80 ,0,2.0));
157  JAKInvMass.at(j).push_back(i.book1D("M12"+tmp,"M_{12,"+tmp+"} (GeV)", 80 ,0,2.0));
158  }
159  }
160  return;
161 }
162 
163 
167  iEvent.getByToken(genparticleCollectionToken_, genParticles );
168 
169  double weight = wmanager_.weight(iEvent);
171  // find taus
172  for (reco::GenParticleCollection::const_iterator iter = genParticles->begin(); iter != genParticles->end(); ++iter) {
173  if(abs(iter->pdgId())==PdtPdgMini::Z0 || abs(iter->pdgId())==PdtPdgMini::Higgs0){
174  spinEffectsZH(&(*iter),weight);
175  }
176  if(abs(iter->pdgId())==15){
177  if(isLastTauinChain(&(*iter))){
178  nTaus->Fill(0.5,weight);
179  int mother = tauMother(&(*iter),weight);
180  if(mother>-1){ // exclude B, D and other non-signal decay modes
181  nPrimeTaus->Fill(0.5,weight);
182  TauPt->Fill(iter->pt(),weight);
183  TauEta->Fill(iter->eta(),weight);
184  TauPhi->Fill(iter->phi(),weight);
185  photons(&(*iter),weight);
187  // Adding JAKID and Mass information
189  unsigned int jak_id, TauBitMask;
190  if(TD.AnalyzeTau(&(*iter),jak_id,TauBitMask,false,false)){
191  JAKID->Fill(jak_id,weight);
192  TauProngs->Fill(TD.nProng(TauBitMask),weight);
193  tauDecayChannel(&(*iter),jak_id,TauBitMask,weight);
194  if(jak_id<=NJAKID){
195  int tcharge=iter->pdgId()/abs(iter->pdgId());
196  std::vector<const reco::GenParticle*> part=TD.Get_TauDecayProducts();
197  spinEffectsWHpm(&(*iter),mother,jak_id,part,weight);
198  TLorentzVector LVQ(0,0,0,0);
199  TLorentzVector LVS12(0,0,0,0);
200  TLorentzVector LVS13(0,0,0,0);
201  TLorentzVector LVS23(0,0,0,0);
202  bool haspart1=false;
203  TVector3 PV,SV;
204  bool hasDL(false);
205  for(unsigned int i=0;i<part.size();i++){
206  if(abs(part.at(i)->pdgId())!=PdtPdgMini::nu_tau && TD.isTauFinalStateParticle(part.at(i)->pdgId()) && !hasDL){
207  TLorentzVector tlv(iter->px(),iter->py(),iter->pz(),iter->energy());
208  PV=TVector3(iter->vx(),iter->vy(),iter->vz());
209  SV=TVector3(part.at(i)->vx(),part.at(i)->vy(),part.at(i)->vz());
210  TVector3 DL=SV-PV;
211  DecayLength->Fill(DL.Dot(tlv.Vect())/tlv.P(),weight);
212  double c(2.99792458E8),Ltau(DL.Mag()/100)/*cm->m*/,beta(iter->p()/iter->mass());
213  LifeTime->Fill( Ltau/(c*beta),weight);
214  }
215 
216  if(TD.isTauFinalStateParticle(part.at(i)->pdgId()) && abs(part.at(i)->pdgId())!=PdtPdgMini::nu_e &&
217  abs(part.at(i)->pdgId())!=PdtPdgMini::nu_mu && abs(part.at(i)->pdgId())!=PdtPdgMini::nu_tau ){
218  TLorentzVector LV(part.at(i)->px(),part.at(i)->py(),part.at(i)->pz(),part.at(i)->energy());
219  LVQ+=LV;
220  if(jak_id==TauDecay::JAK_A1_3PI ||
221  jak_id==TauDecay::JAK_KPIK ||
222  jak_id==TauDecay::JAK_KPIPI
223  ){
224  if((tcharge==part.at(i)->pdgId()/abs(part.at(i)->pdgId()) && TD.nProng(TauBitMask)==3) || (jak_id==TauDecay::JAK_A1_3PI && TD.nProng(TauBitMask)==1 && abs(part.at(i)->pdgId())==PdtPdgMini::pi_plus) ){
225  LVS13+=LV;
226  LVS23+=LV;
227  }
228  else{
229  LVS12+=LV;
230  if(!haspart1 && ((jak_id==TauDecay::JAK_A1_3PI) || (jak_id!=TauDecay::JAK_A1_3PI && abs(part.at(i)->pdgId())==PdtPdgMini::K_plus) )){
231  LVS13+=LV;
232  haspart1=true;
233  }
234  else{
235  LVS23+=LV;
236  }
237  }
238  }
239  }
240  }
241  part.clear();
242  JAKInvMass.at(jak_id).at(0)->Fill(LVQ.M(),weight);
243  if(jak_id==TauDecay::JAK_A1_3PI ||
244  jak_id==TauDecay::JAK_KPIK ||
245  jak_id==TauDecay::JAK_KPIPI
246  ){
247  JAKInvMass.at(jak_id).at(1)->Fill(LVS13.M(),weight);
248  JAKInvMass.at(jak_id).at(2)->Fill(LVS23.M(),weight);
249  JAKInvMass.at(jak_id).at(3)->Fill(LVS12.M(),weight);
250  }
251  }
252  }
253  else{
254  JAKID->Fill(jak_id,weight);
255  }
256  }
257  }
258  }
259  }
260 }//analyze
261 
263  for (unsigned int i=0;i<tau->numberOfMothers();i++) {
264  const reco::GenParticle *mother=static_cast<const reco::GenParticle*>(tau->mother(i));
265  if(mother->pdgId() == tau->pdgId()) return GetMother(mother);
266  return mother;
267  }
268  return tau;
269 }
270 
271 const std::vector<const reco::GenParticle*> TauValidation::GetMothers(const reco::GenParticle* boson){
272  std::vector<const reco::GenParticle*> mothers;
273  for (unsigned int i=0;i<boson->numberOfMothers();i++) {
274  const reco::GenParticle *mother=static_cast<const reco::GenParticle*>(boson->mother(i));
275  if(mother->pdgId() == boson->pdgId()) return GetMothers(mother);
276  mothers.push_back(mother);
277  }
278  return mothers;
279 }
280 
282  return TauValidation::GetMother(tau)->pdgId();
283 }
284 
286  for(unsigned int i = 0; i <tau->numberOfDaughters(); i++){
287  if(tau->daughter(i)->pdgId() == tau->pdgId()) return false;
288  }
289  return true;
290 }
291 
292 void TauValidation::findTauList(const reco::GenParticle* tau,std::vector<const reco::GenParticle*> &TauList){
293  TauList.insert(TauList.begin(),tau);
294  for(unsigned int i=0;i<tau->numberOfMothers();i++) {
295  const reco::GenParticle *mother=static_cast<const reco::GenParticle*>(tau->mother(i));
296  if(mother->pdgId() == tau->pdgId()){
297  findTauList(mother,TauList);
298  }
299  }
300 }
301 
302 void TauValidation::findFSRandBrem(const reco::GenParticle* p, bool doBrem, std::vector<const reco::GenParticle*> &ListofFSR,
303  std::vector<const reco::GenParticle*> &ListofBrem){
304  // 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.
305  if(abs(p->pdgId())==15){
306  if(isLastTauinChain(p)){ doBrem=true;}
307  else{ doBrem=false;}
308  }
309  int photo_ID=22;
310  for(unsigned int i = 0; i <p->numberOfDaughters(); i++){
311  const reco::GenParticle *dau=static_cast<const reco::GenParticle*>(p->daughter(i));
312  if(abs((dau)->pdgId()) == abs(photo_ID) && !doBrem){ListofFSR.push_back(dau);}
313  if(abs((dau)->pdgId()) == abs(photo_ID) && doBrem){ListofBrem.push_back(dau);}
314  if(abs((dau)->pdgId()) != 111 && abs((dau)->pdgId()) != 221){ // remove pi0 and eta decays
315  findFSRandBrem(dau,doBrem,ListofFSR,ListofBrem);
316  }
317  }
318 }
319 
320 void TauValidation::FindPhotosFSR(const reco::GenParticle* p,std::vector<const reco::GenParticle*> &ListofFSR,double &BosonScale){
321  BosonScale=0.0;
322  const reco::GenParticle* m=GetMother(p);
323  int mother_pid=m->pdgId();
324  if(m->pdgId()!=p->pdgId()){
325  for(unsigned int i=0; i <m->numberOfDaughters(); i++){
326  const reco::GenParticle *dau=static_cast<const reco::GenParticle*>(m->daughter(i));
327  if(abs(dau->pdgId()) == 22) {
328  ListofFSR.push_back(dau);
329  }
330  }
331  }
332  if(abs(mother_pid) == 24) BosonScale=1.0; // W
333  if(abs(mother_pid) == 23) BosonScale=2.0; // Z;
334  if(abs(mother_pid) == 22) BosonScale=2.0; // gamma;
335  if(abs(mother_pid) == 25) BosonScale=2.0; // HSM;
336  if(abs(mother_pid) == 35) BosonScale=2.0; // H0;
337  if(abs(mother_pid) == 36) BosonScale=2.0; // A0;
338  if(abs(mother_pid) == 37) BosonScale=1.0; //Hpm;
339 }
340 
342  if(abs(tau->pdgId()) != 15 ) return -3;
343  int mother_pid = findMother(tau);
344  if(mother_pid == -2) return -2;
345  int label = other;
346  if(abs(mother_pid) == 24) label = W;
347  if(abs(mother_pid) == 23) label = Z;
348  if(abs(mother_pid) == 22) label = gamma;
349  if(abs(mother_pid) == 25) label = HSM;
350  if(abs(mother_pid) == 35) label = H0;
351  if(abs(mother_pid) == 36) label = A0;
352  if(abs(mother_pid) == 37) label = Hpm;
353  int mother_shortpid=(abs(mother_pid)%10000);
354  if(mother_shortpid>500 && mother_shortpid<600 )label = B;
355  if(mother_shortpid>400 && mother_shortpid<500)label = D;
356  TauMothers->Fill(label,weight);
357  if(label==B || label == D || label == other) return -1;
358  return mother_pid;
359 }
360 
361 int TauValidation::tauDecayChannel(const reco::GenParticle* tau,int jak_id, unsigned int TauBitMask, double weight){
362  int channel = undetermined;
363  if(tau->status() == 1) channel = stable;
364  int allCount = 0,
365  eCount = 0,
366  muCount = 0,
367  pi0Count = 0,
368  piCount = 0,
369  rhoCount = 0,
370  a1Count = 0,
371  KCount = 0,
372  KstarCount = 0;
373 
374  countParticles(tau,allCount,eCount,muCount,pi0Count,piCount,rhoCount,a1Count,KCount,KstarCount);
375 
376  // resonances
377  if(KCount >= 1) channel = K;
378  if(KstarCount >= 1) channel = Kstar;
379  if(a1Count >= 1) channel = a1;
380  if(rhoCount >= 1) channel = rho;
381  if(channel!=undetermined && weight!=0.0) TauDecayChannels->Fill(channel,weight);
382 
383  // final state products
384  if(piCount == 1 && pi0Count == 0) channel = pi;
385  if(piCount == 1 && pi0Count == 1) channel = pi1pi0;
386  if(piCount == 1 && pi0Count > 1) channel = pinpi0;
387  if(piCount == 3 && pi0Count == 0) channel = tripi;
388  if(piCount == 3 && pi0Count > 0) channel = tripinpi0;
389  if(eCount == 1) channel = electron;
390  if(muCount == 1) channel = muon;
391  if(weight!=0.0) TauDecayChannels->Fill(channel,weight);
392  return channel;
393 }
394 
395 void TauValidation::countParticles(const reco::GenParticle* p,int &allCount, int &eCount, int &muCount,
396  int &pi0Count,int &piCount,int &rhoCount,int &a1Count,int &KCount,int &KstarCount){
397  for(unsigned int i=0; i<p->numberOfDaughters(); i++){
398  const reco::GenParticle *dau=static_cast<const reco::GenParticle*>(p->daughter(i));
399  int pid = dau->pdgId();
400  allCount++;
401  if(abs(pid) == 11) eCount++;
402  if(abs(pid) == 13) muCount++;
403  if(abs(pid) == 111) pi0Count++;
404  if(abs(pid) == 211) piCount++;
405  if(abs(pid) == 213) rhoCount++;
406  if(abs(pid) == 20213) a1Count++;
407  if(abs(pid) == 321) KCount++;
408  if(abs(pid) == 323) KstarCount++;
409  countParticles(dau,allCount,eCount,muCount,pi0Count,piCount,rhoCount,a1Count,KCount,KstarCount);
410  }
411 }
412 
413 
414 
415 void TauValidation::spinEffectsWHpm(const reco::GenParticle* tau,int mother, int decay, std::vector<const reco::GenParticle*> &part,double weight){
416  if(decay == TauDecay::JAK_PION || decay == TauDecay::JAK_MUON || decay == TauDecay::JAK_ELECTRON){ // polarization only for 1-prong hadronic taus with no neutral pions
417  TLorentzVector momP4 = motherP4(tau);
418  TLorentzVector pionP4 = leadingPionP4(tau);
419  pionP4.Boost(-1*momP4.BoostVector());
420  double energy = pionP4.E()/(momP4.M()/2);
421  if(decay == TauDecay::JAK_PION){
422  if(abs(mother) == 24) TauSpinEffectsW_X->Fill(energy,weight);
423  if(abs(mother) == 37) TauSpinEffectsHpm_X->Fill(energy,weight);
424  }
425  if(decay == TauDecay::JAK_MUON){
426  if(abs(mother) == 24) TauSpinEffectsW_muX->Fill(energy,weight);
427  if(abs(mother) == 37) TauSpinEffectsHpm_muX->Fill(energy,weight);
428  }
429  if(decay == TauDecay::JAK_ELECTRON){
430  if(abs(mother) == 24) TauSpinEffectsW_eX->Fill(energy,weight);
431  if(abs(mother) == 37) TauSpinEffectsHpm_eX->Fill(energy,weight);
432  }
433  }
434  else if(decay==TauDecay::JAK_RHO_PIPI0){
435  TLorentzVector rho(0,0,0,0),pi(0,0,0,0);
436  for(unsigned int i=0;i<part.size();i++){
437  TLorentzVector LV(part.at(i)->px(),part.at(i)->py(),part.at(i)->pz(),part.at(i)->energy());
438  if(abs(part.at(i)->pdgId())==PdtPdgMini::pi_plus){pi+=LV; rho+=LV;}
439  if(abs(part.at(i)->pdgId())==PdtPdgMini::pi0){rho+=LV;}
440  }
441  if(abs(mother) == 24) TauSpinEffectsW_UpsilonRho->Fill(2*pi.P()/rho.P()-1,weight);
442  if(abs(mother) == 37) TauSpinEffectsHpm_UpsilonRho->Fill(2*pi.P()/rho.P()-1,weight);
443  }
444  else if(decay==TauDecay::JAK_A1_3PI){ // only for pi2pi0 for now
445  TLorentzVector a1(0,0,0,0),pi_p(0,0,0,0),pi_m(0,0,0,0);
446  int nplus(0),nminus(0);
447  for(unsigned int i=0;i<part.size();i++){
448  TLorentzVector LV(part.at(i)->px(),part.at(i)->py(),part.at(i)->pz(),part.at(i)->energy());
449  if(part.at(i)->pdgId()==PdtPdgMini::pi_plus){ pi_p+=LV; a1+=LV; nplus++;}
450  if(part.at(i)->pdgId()==PdtPdgMini::pi_minus){pi_m+=LV; a1+=LV; nminus++;}
451  }
452  double gamma=0;
453  if(nplus+nminus==3 && nplus==1) gamma=2*pi_p.P()/a1.P()-1;
454  if(nplus+nminus==3 && nminus==1) gamma=2*pi_m.P()/a1.P()-1;
455  else{
456  pi_p+=pi_m; gamma=2*pi_p.P()/a1.P()-1;
457  }
458  if(abs(mother) == 24) TauSpinEffectsW_UpsilonA1->Fill(gamma,weight);
459  if(abs(mother) == 37) TauSpinEffectsHpm_UpsilonA1->Fill(gamma,weight);
460  }
461 }
462 
464  int ntau(0);
465  for(unsigned int i = 0; i <boson->numberOfDaughters(); i++){
466  const reco::GenParticle *dau=static_cast<const reco::GenParticle*>(boson->daughter(i));
467  if(ntau==1 && dau->pdgId() == 15)return;
468  if(boson->pdgId()!= 15 && abs(dau->pdgId()) == 15)ntau++;
469  }
470  if(ntau!=2) return;
471  if(abs(boson->pdgId())==PdtPdgMini::Z0 || abs(boson->pdgId())==PdtPdgMini::Higgs0){
472  TLorentzVector tautau(0,0,0,0);
473  TLorentzVector pipi(0,0,0,0);
474  TLorentzVector taum(0,0,0,0);
475  TLorentzVector taup(0,0,0,0);
476  TLorentzVector rho_plus,rho_minus,pi_rhominus,pi0_rhominus,pi_rhoplus,pi0_rhoplus,pi_plus,pi_minus;
477  bool hasrho_minus(false),hasrho_plus(false),haspi_minus(false),haspi_plus(false);
478  int nSinglePionDecays(0),nSingleMuonDecays(0),nSingleElectronDecays(0);
479  double x1(0),x2(0);
480  TLorentzVector Zboson(boson->px(),boson->py(),boson->pz(),boson->energy());
481  for(unsigned int i = 0; i <boson->numberOfDaughters(); i++){
482  const reco::GenParticle *dau=static_cast<const reco::GenParticle*>(boson->daughter(i));
483  int pid = dau->pdgId();
484  if(abs(findMother(dau)) != 15 && abs(pid) == 15){
486  unsigned int jak_id, TauBitMask;
487  if(TD.AnalyzeTau(dau,jak_id,TauBitMask,false,false)){
488  std::vector<const reco::GenParticle*> part=TD.Get_TauDecayProducts();
489  if(jak_id==TauDecay::JAK_PION || jak_id==TauDecay::JAK_MUON || jak_id==TauDecay::JAK_ELECTRON){
490  if(jak_id==TauDecay::JAK_PION) nSinglePionDecays++;
491  if(jak_id==TauDecay::JAK_MUON) nSingleMuonDecays++;
492  if(jak_id==TauDecay::JAK_ELECTRON) nSingleElectronDecays++;
493  TLorentzVector LVtau(dau->px(),dau->py(),dau->pz(),dau->energy());
494  tautau += LVtau;
495  TLorentzVector LVpi=leadingPionP4(dau);
496  pipi+=LVpi;
497  const HepPDT::ParticleData* pd = fPDGTable->particle(dau->pdgId ());
498  int charge = (int) pd->charge();
499  LVtau.Boost(-1*Zboson.BoostVector());
500  LVpi.Boost(-1*Zboson.BoostVector());
501  if(jak_id==TauDecay::JAK_PION){
502  if(abs(boson->pdgId())==PdtPdgMini::Z0){
503  TauSpinEffectsZ_X->Fill(LVpi.P()/LVtau.E(),weight);
504  if(50.0<Zboson.M() && Zboson.M()<75.0) TauSpinEffectsZ_X50to75->Fill(LVpi.P()/LVtau.E(),weight);
505  if(75.0<Zboson.M() && Zboson.M()<88.0) TauSpinEffectsZ_X75to88->Fill(LVpi.P()/LVtau.E(),weight);
506  if(88.0<Zboson.M() && Zboson.M()<100.0) TauSpinEffectsZ_X88to100->Fill(LVpi.P()/LVtau.E(),weight);
507  if(100.0<Zboson.M() && Zboson.M()<120.0) TauSpinEffectsZ_X100to120->Fill(LVpi.P()/LVtau.E(),weight);
508  if(120.0<Zboson.M()) TauSpinEffectsZ_X120UP->Fill(LVpi.P()/LVtau.E(),weight);
509  }
510  if(abs(boson->pdgId())==PdtPdgMini::Higgs0) TauSpinEffectsH_X->Fill(LVpi.P()/LVtau.E(),weight);
511  }
512  if(charge<0){x1=LVpi.P()/LVtau.E(); taum=LVtau;}
513  else{ x2=LVpi.P()/LVtau.E();}
514  }
515  TLorentzVector LVtau(dau->px(),dau->py(),dau->pz(),dau->energy());
516  if(pid == 15)taum=LVtau;
517  if(pid ==-15)taup=LVtau;
518  if(jak_id==TauDecay::JAK_RHO_PIPI0){
519  for(unsigned int i=0; i<part.size();i++){
520  int pid_d = part.at(i)->pdgId();
521  if(abs(pid_d)==211 || abs(pid_d)==111){
522  TLorentzVector LV(part.at(i)->px(),part.at(i)->py(),part.at(i)->pz(),part.at(i)->energy());
523  if(pid==15){
524  hasrho_minus=true;
525  if(pid_d==-211 ){ pi_rhominus=LV;}
526  if(abs(pid_d)==111 ){ pi0_rhominus=LV;}
527  }
528  if(pid==-15){
529  hasrho_plus=true;
530  if(pid_d==211 ){pi_rhoplus=LV;}
531  if(abs(pid_d)==111 ){pi0_rhoplus=LV;}
532  }
533  }
534  }
535  }
536  if(jak_id==TauDecay::JAK_PION){
537  for(unsigned int i=0; i<part.size();i++){
538  int pid_d = part.at(i)->pdgId();
539  if(abs(pid_d)==211 ){
540  TLorentzVector LV(part.at(i)->px(),part.at(i)->py(),part.at(i)->pz(),part.at(i)->energy());
541  if(pid==15){
542  haspi_minus=true;
543  if(pid_d==-211 ){ pi_minus=LV;}
544  }
545  if(pid==-15){
546  haspi_plus=true;
547  if(pid_d==211 ){pi_plus=LV;}
548  }
549  }
550  }
551  }
552  }
553  }
554  }
555  if(hasrho_minus && hasrho_plus){
556  //compute rhorho
557  rho_minus=pi_rhominus;
558  rho_minus+=pi0_rhominus;
559  rho_plus=pi_rhoplus;
560  rho_plus+=pi0_rhoplus;
561  TLorentzVector rhorho=rho_minus;rhorho+=rho_plus;
562 
563  // boost to rhorho cm
564  TLorentzVector pi_rhoplusb=pi_rhoplus; pi_rhoplusb.Boost(-1*rhorho.BoostVector());
565  TLorentzVector pi0_rhoplusb=pi0_rhoplus; pi0_rhoplusb.Boost(-1*rhorho.BoostVector());
566  TLorentzVector pi_rhominusb=pi_rhominus; pi_rhominusb.Boost(-1*rhorho.BoostVector());
567  TLorentzVector pi0_rhominusb=pi0_rhominus; pi0_rhominusb.Boost(-1*rhorho.BoostVector());
568 
569  // compute n+/-
570  TVector3 n_plus=pi_rhoplusb.Vect().Cross(pi0_rhoplusb.Vect());
571  TVector3 n_minus=pi_rhominusb.Vect().Cross(pi0_rhominusb.Vect());
572 
573  // compute the acoplanarity
574  double Acoplanarity=acos(n_plus.Dot(n_minus)/(n_plus.Mag()*n_minus.Mag()));
575  if(pi_rhominusb.Vect().Dot(n_plus)>0){Acoplanarity*=-1;Acoplanarity+=2*TMath::Pi();}
576 
577  // now boost to tau frame
578  pi_rhoplus.Boost(-1*taup.BoostVector());
579  pi0_rhoplus.Boost(-1*taup.BoostVector());
580  pi_rhominus.Boost(-1*taum.BoostVector());
581  pi0_rhominus.Boost(-1*taum.BoostVector());
582 
583  // compute y1 and y2
584  double y1=(pi_rhoplus.E()-pi0_rhoplus.E())/(pi_rhoplus.E()+pi0_rhoplus.E());
585  double y2=(pi_rhominus.E()-pi0_rhominus.E())/(pi_rhominus.E()+pi0_rhominus.E());
586 
587  // fill histograms
588  if(abs(boson->pdgId())==PdtPdgMini::Higgs0 && y1*y2<0) TauSpinEffectsH_rhorhoAcoplanarityminus->Fill(Acoplanarity,weight);
589  if(abs(boson->pdgId())==PdtPdgMini::Higgs0 && y1*y2>0) TauSpinEffectsH_rhorhoAcoplanarityplus->Fill(Acoplanarity,weight);
590  }
591  if(haspi_minus && haspi_plus){
592  TLorentzVector tauporig=taup;
593  TLorentzVector taumorig=taum;
594 
595  // now boost to Higgs frame
596  pi_plus.Boost(-1*Zboson.BoostVector());
597  pi_minus.Boost(-1*Zboson.BoostVector());
598 
599  taup.Boost(-1*Zboson.BoostVector());
600  taum.Boost(-1*Zboson.BoostVector());
601 
602  if(abs(boson->pdgId())==PdtPdgMini::Higgs0){
603  TauSpinEffectsH_pipiAcollinearity->Fill(acos(pi_plus.Vect().Dot(pi_minus.Vect())/(pi_plus.P()*pi_minus.P())));
604  TauSpinEffectsH_pipiAcollinearityzoom->Fill(acos(pi_plus.Vect().Dot(pi_minus.Vect())/(pi_plus.P()*pi_minus.P())));
605  }
606 
607  double proj_m=taum.Vect().Dot(pi_minus.Vect())/(taum.P()*taum.P());
608  double proj_p=taup.Vect().Dot(pi_plus.Vect())/(taup.P()*taup.P());
609  TVector3 Tau_m=taum.Vect();
610  TVector3 Tau_p=taup.Vect();
611  Tau_m*=proj_m;
612  Tau_p*=proj_p;
613  TVector3 Pit_m=pi_minus.Vect()-Tau_m;
614  TVector3 Pit_p=pi_plus.Vect()-Tau_p;
615 
616  double Acoplanarity=acos(Pit_m.Dot(Pit_p)/(Pit_p.Mag()*Pit_m.Mag()));
617  TVector3 n=Pit_p.Cross(Pit_m);
618  if(n.Dot(Tau_m)/Tau_m.Mag()>0){Acoplanarity*=-1; Acoplanarity+=2*TMath::Pi();}
619  // fill histograms
620  if(abs(boson->pdgId())==PdtPdgMini::Higgs0) TauSpinEffectsH_pipiAcoplanarity->Fill(Acoplanarity,weight);
621  taup=tauporig;
622  taum=taumorig;
623  }
624 
625  if(nSingleMuonDecays==2){
626  if(abs(boson->pdgId())==PdtPdgMini::Z0) TauSpinEffectsZ_muX->Fill(x1,weight);
627  if(abs(boson->pdgId())==PdtPdgMini::Higgs0) TauSpinEffectsH_muX->Fill(x1,weight);
628  }
629  if(nSingleElectronDecays==2){
630  if(abs(boson->pdgId())==PdtPdgMini::Z0) TauSpinEffectsZ_eX->Fill(x1,weight);
631  if(abs(boson->pdgId())==PdtPdgMini::Higgs0) TauSpinEffectsH_eX->Fill(x1,weight);
632  }
633  if(nSinglePionDecays == 2 && tautau.M()!= 0) {
634  for(int i=0;i<zsbins;i++){
635  double zslow=((double)i)*(zsmax-zsmin)/((double)zsbins)+zsmin;
636  double zsup=((double)i+1)*(zsmax-zsmin)/((double)zsbins)+zsmin;
637  double aup=Zstoa(zsup), alow=Zstoa(zslow);
638  if(x2-x1>alow && x2-x1<aup){
639  double zs=(zsup+zslow)/2;
640  if(abs(boson->pdgId())==PdtPdgMini::Z0) TauSpinEffectsZ_Zs->Fill(zs,weight);
641  if(abs(boson->pdgId())==PdtPdgMini::Higgs0) TauSpinEffectsH_Zs->Fill(zs,weight);
642  break;
643  }
644  }
645  if(abs(boson->pdgId())==PdtPdgMini::Z0) TauSpinEffectsZ_MVis->Fill(pipi.M()/tautau.M(),weight);
646  if(abs(boson->pdgId())==PdtPdgMini::Higgs0) TauSpinEffectsH_MVis->Fill(pipi.M()/tautau.M(),weight);
647 
648  if(x1!=0){
649  const std::vector<const reco::GenParticle*> m=GetMothers(boson);
650  int q(0),qbar(0);
651  TLorentzVector Z(0,0,0,0);
652  for(unsigned int i=0;i<m.size();i++){
653  if(m.at(i)->pdgId()==PdtPdgMini::d || m.at(i)->pdgId()==PdtPdgMini::u ){q++;}
654  if(m.at(i)->pdgId()==PdtPdgMini::anti_d || m.at(i)->pdgId()==PdtPdgMini::anti_u ){qbar++;}
655  }
656  if(q==1 && qbar==1){// assume q has largest E (valence vs see quarks)
657  if(taum.Vect().Dot(Zboson.Vect())/(Zboson.P()*taum.P())>0){
658  if(abs(boson->pdgId())==PdtPdgMini::Z0) TauSpinEffectsZ_Xf->Fill(x1,weight);
659  if(abs(boson->pdgId())==PdtPdgMini::Higgs0) TauSpinEffectsH_Xf->Fill(x1,weight);
660  }
661  else{
662  if(abs(boson->pdgId())==PdtPdgMini::Z0) TauSpinEffectsZ_Xb->Fill(x1,weight);
663  if(abs(boson->pdgId())==PdtPdgMini::Higgs0) TauSpinEffectsH_Xb->Fill(x1,weight);
664  }
665  }
666  }
667  }
668  }
669 }
670 
671 double TauValidation::Zstoa(double zs){
672  double a=1-sqrt(fabs(1.0-2*fabs(zs)));
673  if(zs<0){
674  a*=-1.0;
675  }
676  return a;
677 }
678 
679 
681  return leadingPionP4(tau).P();
682 }
683 
685  TLorentzVector p4(0,0,0,0);
686  for(unsigned int i = 0; i <tau->numberOfDaughters(); i++){
687  const reco::GenParticle *dau=static_cast<const reco::GenParticle*>(tau->daughter(i));
688  int pid = dau->pdgId();
689  if(abs(pid) == 15) return leadingPionP4(dau);
690  if(!(abs(pid)==211 || abs(pid)==13 || abs(pid)==11)) continue;
691  if(dau->p() > p4.P()) p4 = TLorentzVector(dau->px(),dau->py(),dau->pz(),dau->energy());
692  }
693  return p4;
694 }
695 
697  const reco::GenParticle* m=GetMother(tau);
698  return TLorentzVector(m->px(),m->py(),m->pz(),m->energy());
699 }
700 
702  TLorentzVector p4(tau->px(),tau->py(),tau->pz(),tau->energy());
703  for(unsigned int i = 0; i <tau->numberOfDaughters(); i++){
704  const reco::GenParticle *dau=static_cast<const reco::GenParticle*>(tau->daughter(i));
705  int pid = dau->pdgId();
706  if(abs(pid) == 15) return visibleTauEnergy(dau);
707  if(abs(pid) == 12 || abs(pid) == 14 || abs(pid) == 16) {
708  p4-=TLorentzVector(dau->px(),dau->py(),dau->pz(),dau->energy());
709  }
710  }
711  return p4.E();
712 }
713 
715  // Find First tau in chain
716  std::vector<const reco::GenParticle*> TauList;
717  findTauList(tau,TauList);
718 
719  // Get List of Gauge Boson to tau(s) FSR and Brem
720  bool passedW=false;
721  std::vector<const reco::GenParticle*> ListofFSR; ListofFSR.clear();
722  std::vector<const reco::GenParticle*> ListofBrem; ListofBrem.clear();
723  std::vector<const reco::GenParticle*> FSR_photos; FSR_photos.clear();
724  double BosonScale(1);
725  if(TauList.size()>0){
726  TauValidation::findFSRandBrem(TauList.at(0),passedW,ListofFSR,ListofBrem);
727  TauValidation::FindPhotosFSR(TauList.at(0),FSR_photos,BosonScale);
728 
729  // Add the Tau Brem. information
730  TauBremPhotonsN->Fill(ListofBrem.size(),weight);
731  double photonPtSum=0;
732  for(unsigned int i=0;i<ListofBrem.size();i++){
733  photonPtSum+=ListofBrem.at(i)->pt();
734  TauBremPhotonsPt->Fill(ListofBrem.at(i)->pt(),weight);
735  }
736  TauBremPhotonsPtSum->Fill(photonPtSum,weight);
737 
738  // Now add the Gauge Boson FSR information
739  if(BosonScale!=0){
740  TauFSRPhotonsN->Fill(ListofFSR.size(),weight);
741  photonPtSum=0;
742  for(unsigned int i=0;i<ListofFSR.size();i++){
743  photonPtSum+=ListofFSR.at(i)->pt();
744  TauFSRPhotonsPt->Fill(ListofFSR.at(i)->pt(),weight);
745  }
746  double FSR_photosSum(0);
747  for(unsigned int i=0;i<FSR_photos.size();i++){
748  FSR_photosSum+=FSR_photos.at(i)->pt();
749  TauFSRPhotonsPt->Fill(FSR_photos.at(i)->pt()/BosonScale,weight*BosonScale);
750  }
751  TauFSRPhotonsPtSum->Fill(photonPtSum+FSR_photosSum/BosonScale,weight);
752  }
753  }
754 }
755 
MonitorElement * TauFSRPhotonsN
Definition: TauValidation.h:95
const double beta
const double Pi
virtual double energy() const GCC11_FINAL
energy
int i
Definition: DBlmapReader.cc:9
int findMother(const reco::GenParticle *)
virtual ~TauValidation()
MonitorElement * TauSpinEffectsH_rhorhoAcoplanarityplus
Definition: TauValidation.h:95
MonitorElement * TauPhi
Definition: TauValidation.h:95
MonitorElement * TauSpinEffectsW_X
Definition: TauValidation.h:95
MonitorElement * TauSpinEffectsHpm_UpsilonA1
Definition: TauValidation.h:95
std::vector< std::vector< MonitorElement * > > JAKInvMass
MonitorElement * TauSpinEffectsZ_Xb
Definition: TauValidation.h:95
MonitorElement * nTaus
Definition: TauValidation.h:94
virtual double p() const GCC11_FINAL
magnitude of momentum vector
WeightManager wmanager_
Definition: TauValidation.h:64
MonitorElement * DecayLength
Definition: TauValidation.h:95
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
virtual int pdgId() const GCC11_FINAL
PDG identifier.
TauValidation(const edm::ParameterSet &)
MonitorElement * TauPt
Definition: TauValidation.h:95
MonitorElement * TauSpinEffectsH_MVis
Definition: TauValidation.h:95
MonitorElement * TauMothers
Definition: TauValidation.h:95
MonitorElement * TauSpinEffectsHpm_eX
Definition: TauValidation.h:95
MonitorElement * TauSpinEffectsZ_X100to120
Definition: TauValidation.h:95
void FindPhotosFSR(const reco::GenParticle *p, std::vector< const reco::GenParticle * > &ListofFSR, double &BosonScale)
MonitorElement * TauSpinEffectsW_muX
Definition: TauValidation.h:95
MonitorElement * TauBremPhotonsN
Definition: TauValidation.h:95
void spinEffectsWHpm(const reco::GenParticle *, int, int, std::vector< const reco::GenParticle * > &part, double weight)
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
Definition: TauValidation.h:95
bool isLastTauinChain(const reco::GenParticle *tau)
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
MonitorElement * TauBremPhotonsPtSum
Definition: TauValidation.h:95
MonitorElement * TauSpinEffectsH_pipiAcoplanarity
Definition: TauValidation.h:95
MonitorElement * TauSpinEffectsHpm_UpsilonRho
Definition: TauValidation.h:95
virtual double pz() const GCC11_FINAL
z coordinate of momentum vector
virtual double py() const GCC11_FINAL
y coordinate of momentum vector
double leadingPionMomentum(const reco::GenParticle *, double weight)
void getData(T &iHolder) const
Definition: EventSetup.h:78
MonitorElement * TauSpinEffectsZ_eX
Definition: TauValidation.h:95
void Fill(long long x)
MonitorElement * TauProngs
Definition: TauValidation.h:95
TLorentzVector motherP4(const reco::GenParticle *)
bool isTauFinalStateParticle(int pdgid)
Definition: TauDecay.cc:36
virtual void dqmBeginRun(const edm::Run &r, const edm::EventSetup &c) override
MonitorElement * TauSpinEffectsH_pipiAcollinearityzoom
Definition: TauValidation.h:95
math::XYZTLorentzVectorD LV
void photons(const reco::GenParticle *, double weight)
int iEvent
Definition: GenABIO.cc:230
MonitorElement * TauSpinEffectsW_UpsilonRho
Definition: TauValidation.h:95
const reco::GenParticle * GetMother(const reco::GenParticle *tau)
MonitorElement * TauSpinEffectsZ_X88to100
Definition: TauValidation.h:95
MonitorElement * TauSpinEffectsH_rhorhoAcoplanarityminus
Definition: TauValidation.h:95
virtual int status() const GCC11_FINAL
status word
virtual void analyze(edm::Event const &, edm::EventSetup const &) override
virtual size_t numberOfMothers() const
number of mothers
MonitorElement * TauSpinEffectsZ_X50to75
Definition: TauValidation.h:95
T sqrt(T t)
Definition: SSEVec.h:48
double p4[4]
Definition: TauolaWrapper.h:92
virtual size_t numberOfDaughters() const
number of daughters
virtual const Candidate * daughter(size_type) const
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
MonitorElement * TauSpinEffectsH_Xf
Definition: TauValidation.h:95
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
MonitorElement * TauFSRPhotonsPt
Definition: TauValidation.h:95
virtual double px() const GCC11_FINAL
x coordinate of momentum vector
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
MonitorElement * TauSpinEffectsZ_muX
Definition: TauValidation.h:95
MonitorElement * TauDecayChannels
Definition: TauValidation.h:95
edm::InputTag hepmcCollection_
Definition: TauValidation.h:89
MonitorElement * TauSpinEffectsH_Zs
Definition: TauValidation.h:95
MonitorElement * TauFSRPhotonsPtSum
Definition: TauValidation.h:95
HepPDT::ParticleData ParticleData
int tauMother(const reco::GenParticle *, double weight)
virtual void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
std::vector< const reco::GenParticle * > Get_TauDecayProducts()
const std::vector< const reco::GenParticle * > GetMothers(const reco::GenParticle *boson)
MonitorElement * TauSpinEffectsZ_MVis
Definition: TauValidation.h:95
virtual int pdgId() const =0
PDG identifier.
MonitorElement * TauSpinEffectsZ_X75to88
Definition: TauValidation.h:95
double visibleTauEnergy(const reco::GenParticle *)
int tauDecayChannel(const reco::GenParticle *tau, int jak_id, unsigned int TauBitMask, double weight)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * TauSpinEffectsH_pipiAcollinearity
Definition: TauValidation.h:95
part
Definition: HCALResponse.h:20
tuple pid
Definition: sysUtil.py:22
MonitorElement * TauSpinEffectsZ_X
Definition: TauValidation.h:95
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.
Definition: TauValidation.h:92
bool AnalyzeTau(const reco::GenParticle *Tau, unsigned int &JAK_ID, unsigned int &TauBitMask, bool dores, bool dopi0)
MonitorElement * TauSpinEffectsH_eX
Definition: TauValidation.h:95
MonitorElement * TauSpinEffectsW_eX
Definition: TauValidation.h:95
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
edm::EDGetTokenT< reco::GenParticleCollection > genparticleCollectionToken_
double Zstoa(double zs)
void findTauList(const reco::GenParticle *tau, std::vector< const reco::GenParticle * > &TauList)
void findFSRandBrem(const reco::GenParticle *p, bool doBrem, std::vector< const reco::GenParticle * > &ListofFSR, std::vector< const reco::GenParticle * > &ListofBrem)
MonitorElement * LifeTime
Definition: TauValidation.h:95
MonitorElement * TauSpinEffectsH_muX
Definition: TauValidation.h:95
MonitorElement * TauSpinEffectsW_UpsilonA1
Definition: TauValidation.h:95
double a
Definition: hdecay.h:121
void countParticles(const reco::GenParticle *p, int &allCount, int &eCount, int &muCount, int &pi0Count, int &piCount, int &rhoCount, int &a1Count, int &KCount, int &KstarCount)
MonitorElement * TauEta
Definition: TauValidation.h:95
edm::InputTag genparticleCollection_
Definition: TauValidation.h:88
unsigned int NJAKID
MonitorElement * TauSpinEffectsZ_Zs
Definition: TauValidation.h:95
MonitorElement * TauSpinEffectsHpm_X
Definition: TauValidation.h:95
int weight
Definition: histoStyle.py:50
MonitorElement * TauSpinEffectsHpm_muX
Definition: TauValidation.h:95
double weight(const edm::Event &)
MonitorElement * TauSpinEffectsZ_X120UP
Definition: TauValidation.h:95
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
MonitorElement * TauSpinEffectsH_Xb
Definition: TauValidation.h:95
virtual const Candidate * mother(size_type=0) const
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
MonitorElement * nPrimeTaus
Definition: TauValidation.h:94
MonitorElement * TauBremPhotonsPt
Definition: TauValidation.h:95
void spinEffectsZH(const reco::GenParticle *boson, double weight)
Definition: Run.h:41
unsigned int nProng(unsigned int &TauBitMask)
Definition: TauDecay.cc:256
MonitorElement * TauSpinEffectsH_X
Definition: TauValidation.h:95
TLorentzVector leadingPionP4(const reco::GenParticle *)