#include <DTTrigTest.h>
Public Member Functions | |
void | analyze (const edm::Event &iEvent, const edm::EventSetup &iEventSetup) |
Analyze function executed on all the events. | |
void | beginJob () |
Create tree and Branches. | |
void | beginRun (const edm::Run &iRun, const edm::EventSetup &iEventSetup) |
Create DTTrig instance and TUs. | |
DTTrigTest (const edm::ParameterSet &pset) | |
Constructor. | |
void | endJob () |
Close Tree and write File. | |
~DTTrigTest () | |
Destructor. | |
Private Attributes | |
int | bbx [100] |
int | bcod [100] |
float | bdirx [100] |
float | bdiry [100] |
float | bdirz [100] |
int | bk [100] |
int | bnum [100] |
float | bposx [100] |
float | bposy [100] |
float | bposz [100] |
int | bsect [100] |
int | bsl [100] |
int | bstat [100] |
int | bwh [100] |
int | bx [100] |
int | chagen [10] |
float | etagen [10] |
int | eventn |
bool | my_debug |
edm::ParameterSet | my_params |
TFile * | my_rootfile |
TTree * | my_tree |
DTTrig * | my_trig |
int | nbti |
int | ngen |
int | nscphi |
int | nsctheta |
int | ntraco |
int | ntsphi |
int | ntstheta |
float | phigen [10] |
float | ptgen [10] |
float | pxgen [10] |
float | pygen [10] |
float | pzgen [10] |
int | runn |
int | sbx [40] |
int | scod [40] |
int | scphbx [40] |
int | scphcod [40] |
float | scphdirx [100] |
float | scphdiry [100] |
float | scphdirz [100] |
int | scphphi [40] |
int | scphphib [40] |
float | scphposx [100] |
float | scphposy [100] |
float | scphposz [100] |
int | scphsect [40] |
int | scphstat [40] |
int | scphwh [40] |
int | scthbx [40] |
int | scthcode [40][7] |
int | scthpos [40][7] |
int | scthqual [40][7] |
int | scthsect [40] |
int | scthstat [40] |
int | scthwh [40] |
float | sdirx [100] |
float | sdiry [100] |
float | sdirz [100] |
int | sphi [40] |
int | sphib [40] |
float | sposx [100] |
float | sposy [100] |
float | sposz [100] |
int | ssect [40] |
int | sstat [40] |
int | swh [40] |
int | tbx [80] |
int | tcod [80] |
float | tdirx [100] |
float | tdiry [100] |
float | tdirz [100] |
int | thbx [40] |
int | thcode [40][7] |
int | thpos [40][7] |
int | thqual [40][7] |
int | thsect [40] |
int | thstat [40] |
int | thwh [40] |
int | tk [80] |
int | tnum [80] |
float | tposx [100] |
float | tposy [100] |
float | tposz [100] |
int | tsect [80] |
int | tstat [80] |
int | twh [80] |
int | tx [80] |
float | vxgen [10] |
float | vygen [10] |
float | vzgen [10] |
float | weight |
Static Private Attributes | |
static const double | my_TtoTDC = 32./25. |
EDAnalyzer that generates a rootfile useful for L1-DTTrigger debugging and performance studies
EDAnalyzer that generates a rootfile useful for L1-DTTrigger debugging and performance studies
Definition at line 34 of file DTTrigTest.h.
DTTrigTest::DTTrigTest | ( | const edm::ParameterSet & | pset | ) |
Constructor.
Definition at line 53 of file DTTrigTest.cc.
References gather_cfg::cout, edm::ParameterSet::getUntrackedParameter(), my_debug, my_params, my_rootfile, my_tree, and estimatePileup_makeJSON::outputfile.
: my_trig(0) { my_debug= pset.getUntrackedParameter<bool>("debug"); string outputfile = pset.getUntrackedParameter<string>("outputFileName"); if (my_debug) cout << "[DTTrigTest] Creating rootfile " << outputfile <<endl; my_rootfile = new TFile(outputfile.c_str(),"RECREATE"); my_tree = new TTree("h1","GMT",0); my_params = pset; if (my_debug) cout << "[DTTrigTest] Constructor executed!!!" << endl; }
DTTrigTest::~DTTrigTest | ( | ) |
Destructor.
Definition at line 67 of file DTTrigTest.cc.
References gather_cfg::cout, my_debug, my_rootfile, and my_trig.
{ if (my_trig != 0) delete my_trig; delete my_rootfile; if (my_debug) cout << "[DTTrigTest] Destructor executed!!!" << endl; }
void DTTrigTest::analyze | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iEventSetup | ||
) | [virtual] |
Analyze function executed on all the events.
Implements edm::EDAnalyzer.
Definition at line 216 of file DTTrigTest.cc.
References abs, bbx, bcod, bdirx, bdiry, bdirz, bk, bnum, bposx, bposy, bposz, bsect, bsl, bstat, DTTrig::BtiTrigs(), bwh, bx, chagen, DeDxDiscriminatorTools::charge(), DTTrig::CMSDirection(), DTTrig::CMSPosition(), gather_cfg::cout, dir, eta(), etagen, edm::EventID::event(), eventn, edm::Event::getByLabel(), i, edm::EventBase::id(), M_PI, MAXGEN, my_debug, my_tree, my_trig, nbti, ngen, nscphi, nsctheta, ntraco, ntsphi, ntstheta, phi, phigen, pos, ptgen, pxgen, pygen, pzgen, edm::EventID::run(), runn, sbx, scod, scphbx, scphcod, scphdirx, scphdiry, scphdirz, scphphi, scphphib, scphposx, scphposy, scphposz, scphsect, scphstat, DTTrig::SCPhTrigs(), scphwh, scthbx, scthcode, scthpos, scthqual, scthsect, scthstat, DTTrig::SCThTrigs(), scthwh, sdirx, sdiry, sdirz, sphi, sphib, sposx, sposy, sposz, ssect, sstat, swh, tbx, tcod, tdirx, tdiry, tdirz, thbx, thcode, thpos, thqual, thsect, thstat, thwh, tk, tnum, tposx, tposy, tposz, DTTrig::TracoTrigs(), DTTrig::triggerReco(), tsect, DTTrig::TSPhTrigs(), tstat, DTTrig::TSThTrigs(), twh, tx, vxgen, vygen, vzgen, weight, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().
{ const int MAXGEN = 10; const float ptcut = 1.0; const float etacut = 2.4; my_trig->triggerReco(iEvent,iEventSetup); if (my_debug) cout << "[DTTrigTest] Trigger algorithm executed for run " << iEvent.id().run() <<" event " << iEvent.id().event() << endl; // GENERAL Block runn = iEvent.id().run(); eventn = iEvent.id().event(); weight = 1; // FIXME what to do with this varable? // GEANT Block Handle<vector<SimTrack> > MyTracks; Handle<vector<SimVertex> > MyVertexes; iEvent.getByLabel("g4SimHits",MyTracks); iEvent.getByLabel("g4SimHits",MyVertexes); vector<SimTrack>::const_iterator itrack; ngen=0; if (my_debug) cout << "[DTTrigTest] Tracks found in the detector (not only muons) " << MyTracks->size() <<endl; for (itrack=MyTracks->begin(); itrack!=MyTracks->end(); itrack++){ if ( abs(itrack->type())==13){ math::XYZTLorentzVectorD momentum = itrack->momentum(); float pt = momentum.Pt(); float eta = momentum.eta(); if ( pt>ptcut && fabs(eta)<etacut ){ float phi = momentum.phi(); int charge = static_cast<int> (-itrack->type()/13); //static_cast<int> (itrack->charge()); if ( phi<0 ) phi = 2*M_PI + phi; int vtxindex = itrack->vertIndex(); float gvx=0,gvy=0,gvz=0; if (vtxindex >-1){ gvx=MyVertexes->at(vtxindex).position().x(); gvy=MyVertexes->at(vtxindex).position().y(); gvz=MyVertexes->at(vtxindex).position().z(); } if ( ngen < MAXGEN ) { pxgen[ngen]=momentum.x(); pygen[ngen]=momentum.y(); pzgen[ngen]=momentum.z(); ptgen[ngen]=pt; etagen[ngen]=eta; phigen[ngen]=phi; chagen[ngen]=charge; vxgen[ngen]=gvx; vygen[ngen]=gvy; vzgen[ngen]=gvz; ngen++; } } } } // L1 Local Trigger Block // BTI vector<DTBtiTrigData> btitrigs = my_trig->BtiTrigs(); vector<DTBtiTrigData>::const_iterator pbti; int ibti = 0; if (my_debug) cout << "[DTTrigTest] " << btitrigs.size() << " BTI triggers found" << endl; for ( pbti = btitrigs.begin(); pbti != btitrigs.end(); pbti++ ) { if ( ibti < 100 ) { bwh[ibti]=pbti->wheel(); bstat[ibti]=pbti->station(); bsect[ibti]=pbti->sector(); bsl[ibti]=pbti->btiSL(); bnum[ibti]=pbti->btiNumber(); bbx[ibti]=pbti->step(); bcod[ibti]=pbti->code(); bk[ibti]=pbti->K(); bx[ibti]=pbti->X(); GlobalPoint pos = my_trig->CMSPosition(&(*pbti)); GlobalVector dir = my_trig->CMSDirection(&(*pbti)); bposx[ibti] = pos.x(); bposy[ibti] = pos.y(); bposz[ibti] = pos.z(); bdirx[ibti] = dir.x(); bdiry[ibti] = dir.y(); bdirz[ibti] = dir.z(); ibti++; } } nbti = ibti; //TRACO vector<DTTracoTrigData> tracotrigs = my_trig->TracoTrigs(); vector<DTTracoTrigData>::const_iterator ptc; int itraco = 0; if (my_debug) cout << "[DTTrigTest] " << tracotrigs.size() << " TRACO triggers found" << endl; for (ptc=tracotrigs.begin(); ptc!=tracotrigs.end(); ptc++) { if (itraco<80) { twh[itraco]=ptc->wheel(); tstat[itraco]=ptc->station(); tsect[itraco]=ptc->sector(); tnum[itraco]=ptc->tracoNumber(); tbx[itraco]=ptc->step(); tcod[itraco]=ptc->code(); tk[itraco]=ptc->K(); tx[itraco]=ptc->X(); GlobalPoint pos = my_trig->CMSPosition(&(*ptc)); GlobalVector dir = my_trig->CMSDirection(&(*ptc)); tposx[itraco] = pos.x(); tposy[itraco] = pos.y(); tposz[itraco] = pos.z(); tdirx[itraco] = dir.x(); tdiry[itraco] = dir.y(); tdirz[itraco] = dir.z(); itraco++; } } ntraco = itraco; //TSPHI vector<DTChambPhSegm> tsphtrigs = my_trig->TSPhTrigs(); vector<DTChambPhSegm>::const_iterator ptsph; int itsphi = 0; if (my_debug ) cout << "[DTTrigTest] " << tsphtrigs.size() << " TSPhi triggers found" << endl; for (ptsph=tsphtrigs.begin(); ptsph!=tsphtrigs.end(); ptsph++) { if (itsphi<40 ) { swh[itsphi] = ptsph->wheel(); sstat[itsphi] = ptsph->station(); ssect[itsphi] = ptsph->sector(); sbx[itsphi] = ptsph->step(); scod[itsphi] = ptsph->oldCode(); sphi[itsphi] = ptsph->phi(); sphib[itsphi] = ptsph->phiB(); GlobalPoint pos = my_trig->CMSPosition(&(*ptsph)); GlobalVector dir = my_trig->CMSDirection(&(*ptsph)); sposx[itsphi] = pos.x(); sposy[itsphi] = pos.y(); sposz[itsphi] = pos.z(); sdirx[itsphi] = dir.x(); sdiry[itsphi] = dir.y(); sdirz[itsphi] = dir.z(); itsphi++; } } ntsphi = itsphi; //TSTHETA vector<DTChambThSegm> tsthtrigs = my_trig->TSThTrigs(); vector<DTChambThSegm>::const_iterator ptsth; int itstheta = 0; if (my_debug) cout << "[DTTrigTest] " << tsthtrigs.size() << " TSTheta triggers found" << endl; for (ptsth=tsthtrigs.begin(); ptsth!=tsthtrigs.end(); ptsth++) { if (itstheta<40 ) { thwh[itstheta] = ptsth->ChamberId().wheel(); thstat[itstheta] = ptsth->ChamberId().station(); thsect[itstheta] = ptsth->ChamberId().sector(); thbx[itstheta] = ptsth->step(); for(int i=0;i<7;i++) { thcode[itstheta][i] = ptsth->code(i); thpos[itstheta][i] = ptsth->position(i); thqual[itstheta][i] = ptsth->quality(i); } itstheta++; } } ntstheta = itstheta; //SCPHI vector<DTSectCollPhSegm> scphtrigs = my_trig->SCPhTrigs(); vector<DTSectCollPhSegm>::const_iterator pscph; int iscphi = 0; if (my_debug ) cout << "[DTTrigTest] " << scphtrigs.size() << " SectCollPhi triggers found" << endl; for (pscph=scphtrigs.begin(); pscph!=scphtrigs.end(); pscph++) { if (iscphi<40 ) { const DTChambPhSegm *seg = (*pscph).tsPhiTrig(); scphwh[iscphi] = pscph->wheel(); scphstat[iscphi] = pscph->station(); scphsect[iscphi] = pscph->sector(); scphbx[iscphi] = pscph->step(); scphcod[iscphi] = pscph->oldCode(); scphphi[iscphi] = pscph->phi(); scphphib[iscphi] = pscph->phiB(); GlobalPoint pos = my_trig->CMSPosition(seg); GlobalVector dir = my_trig->CMSDirection(seg); scphposx[iscphi] = pos.x(); scphposy[iscphi] = pos.y(); scphposz[iscphi] = pos.z(); scphdirx[iscphi] = dir.x(); scphdiry[iscphi] = dir.y(); scphdirz[iscphi] = dir.z(); iscphi++; } } nscphi = iscphi; //SCTHETA vector<DTSectCollThSegm> scthtrigs = my_trig->SCThTrigs(); vector<DTSectCollThSegm>::const_iterator pscth; int isctheta = 0; if (my_debug) cout << "[DTTrigTest] " << scthtrigs.size() << " SectCollTheta triggers found" << endl; for (pscth=scthtrigs.begin(); pscth!=scthtrigs.end(); pscth++) { if (isctheta<40 ) { scthwh[isctheta] = pscth->ChamberId().wheel(); scthstat[isctheta] = pscth->ChamberId().station(); scthsect[isctheta] = pscth->ChamberId().sector(); scthbx[isctheta] = pscth->step(); for(int i=0;i<7;i++) { scthcode[isctheta][i] = pscth->code(i); scthpos[isctheta][i] = pscth->position(i); scthqual[isctheta][i] = pscth->quality(i); } isctheta++; } } nsctheta = isctheta; //Fill the tree my_tree->Fill(); }
void DTTrigTest::beginJob | ( | void | ) | [virtual] |
Create tree and Branches.
Reimplemented from edm::EDAnalyzer.
Definition at line 87 of file DTTrigTest.cc.
References bbx, bcod, bdirx, bdiry, bdirz, bk, bnum, bposx, bposy, bposz, bsect, bsl, bstat, bwh, bx, chagen, etagen, eventn, my_tree, nbti, ngen, nscphi, nsctheta, ntraco, ntsphi, ntstheta, phigen, ptgen, pxgen, pygen, pzgen, runn, sbx, scod, scphbx, scphcod, scphdirx, scphdiry, scphdirz, scphphi, scphphib, scphposx, scphposy, scphposz, scphsect, scphstat, scphwh, scthbx, scthcode, scthpos, scthqual, scthsect, scthstat, scthwh, sdirx, sdiry, sdirz, sphi, sphib, sposx, sposy, sposz, ssect, sstat, swh, tbx, tcod, tdirx, tdiry, tdirz, thbx, thcode, thpos, thqual, thsect, thstat, thwh, tk, tnum, tposx, tposy, tposz, tsect, tstat, twh, tx, vxgen, vygen, vzgen, and weight.
{ // get DTConfigManager // ESHandle< DTConfigManager > confManager ; // iEventSetup.get< DTConfigManagerRcd >().get( confManager ) ; //for testing purpose.... //DTBtiId btiid(1,1,1,1,1); //confManager->getDTConfigBti(btiid)->print(); // my_trig = new DTTrig(my_params); // my_trig->createTUs(iEventSetup); // if (my_debug) // cout << "[DTTrigTest] TU's Created" << endl; // BOOKING of the tree's varables // GENERAL block branches my_tree->Branch("Run",&runn,"Run/I"); my_tree->Branch("Event",&eventn,"Event/I"); my_tree->Branch("Weight",&weight,"Weight/F"); // GEANT block branches my_tree->Branch("Ngen",&ngen,"Ngen/I"); my_tree->Branch("Pxgen",pxgen,"Pxgen[Ngen]/F"); my_tree->Branch("Pygen",pygen,"Pygen[Ngen]/F"); my_tree->Branch("Pzgen",pzgen,"Pzgen[Ngen]/F"); my_tree->Branch("Ptgen",ptgen,"Ptgen[Ngen]/F"); my_tree->Branch("Etagen",etagen,"Etagen[Ngen]/F"); my_tree->Branch("Phigen",phigen,"Phigen[Ngen]/F"); my_tree->Branch("Chagen",chagen,"Chagen[Ngen]/I"); my_tree->Branch("Vxgen",vxgen,"Vxgen[Ngen]/F"); my_tree->Branch("Vygen",vygen,"Vygen[Ngen]/F"); my_tree->Branch("Vzgen",vzgen,"Vzgen[Ngen]/F"); // L1MuDTBtiChipS block my_tree->Branch("Nbti",&nbti,"Nbti/I"); my_tree->Branch("bwh",bwh,"bwh[Nbti]/I"); my_tree->Branch("bstat",bstat,"bstat[Nbti]/I"); my_tree->Branch("bsect",bsect,"bsect[Nbti]/I"); my_tree->Branch("bsl",bsl,"bsl[Nbti]/I"); my_tree->Branch("bnum",bnum,"bnum[Nbti]/I"); my_tree->Branch("bbx",bbx,"bbx[Nbti]/I"); my_tree->Branch("bcod",bcod,"bcod[Nbti]/I"); my_tree->Branch("bk",bk,"bk[Nbti]/I"); my_tree->Branch("bx",bx,"bx[Nbti]/I"); my_tree->Branch("bposx",bposx,"bposx[Nbti]/F"); my_tree->Branch("bposy",bposy,"bposy[Nbti]/F"); my_tree->Branch("bposz",bposz,"bposz[Nbti]/F"); my_tree->Branch("bdirx",bdirx,"bdirx[Nbti]/F"); my_tree->Branch("bdiry",bdiry,"bdiry[Nbti]/F"); my_tree->Branch("bdirz",bdirz,"bdirz[Nbti]/F"); // L1MuDTTracoChipS block my_tree->Branch("Ntraco",&ntraco,"Ntraco/I"); my_tree->Branch("twh",twh,"twh[Ntraco]/I"); my_tree->Branch("tstat",tstat,"tstat[Ntraco]/I"); my_tree->Branch("tsect",tsect,"tsect[Ntraco]/I"); my_tree->Branch("tnum",tnum,"tnum[Ntraco]/I"); my_tree->Branch("tbx",tbx,"tbx[Ntraco]/I"); my_tree->Branch("tcod",tcod,"tcod[Ntraco]/I"); my_tree->Branch("tk",tk,"tk[Ntraco]/I"); my_tree->Branch("tx",tx,"tx[Ntraco]/I"); my_tree->Branch("tposx",tposx,"tposx[Ntraco]/F"); my_tree->Branch("tposy",tposy,"tposy[Ntraco]/F"); my_tree->Branch("tposz",tposz,"tposz[Ntraco]/F"); my_tree->Branch("tdirx",tdirx,"tdirx[Ntraco]/F"); my_tree->Branch("tdiry",tdiry,"tdiry[Ntraco]/F"); my_tree->Branch("tdirz",tdirz,"tdirz[Ntraco]/F"); // TSPHI block my_tree->Branch("Ntsphi",&ntsphi,"Ntsphi/I"); my_tree->Branch("swh",swh,"swh[Ntsphi]/I"); my_tree->Branch("sstat",sstat,"sstat[Ntsphi]/I"); my_tree->Branch("ssect",ssect,"ssect[Ntsphi]/I"); my_tree->Branch("sbx",sbx,"sbx[Ntsphi]/I"); my_tree->Branch("scod",scod,"scod[Ntsphi]/I"); my_tree->Branch("sphi",sphi,"sphi[Ntsphi]/I"); my_tree->Branch("sphib",sphib,"sphib[Ntsphi]/I"); my_tree->Branch("sposx",sposx,"sposx[Ntsphi]/F"); my_tree->Branch("sposy",sposy,"sposy[Ntsphi]/F"); my_tree->Branch("sposz",sposz,"sposz[Ntsphi]/F"); my_tree->Branch("sdirx",sdirx,"sdirx[Ntsphi]/F"); my_tree->Branch("sdiry",sdiry,"sdiry[Ntsphi]/F"); my_tree->Branch("sdirz",sdirz,"sdirz[Ntsphi]/F"); // TSTHETA block my_tree->Branch("Ntstheta",&ntstheta,"Ntstheta/I"); my_tree->Branch("thwh",thwh,"thwh[Ntstheta]/I"); my_tree->Branch("thstat",thstat,"thstat[Ntstheta]/I"); my_tree->Branch("thsect",thsect,"thsect[Ntstheta]/I"); my_tree->Branch("thbx",thbx,"thbx[Ntstheta]/I"); my_tree->Branch("thcode",thcode,"thcode[Ntstheta][7]/I"); my_tree->Branch("thpos",thpos,"thpos[Ntstheta][7]/I"); my_tree->Branch("thqual",thqual,"thqual[Ntstheta][7]/I"); // SC PHI block my_tree->Branch("Nscphi",&nscphi,"Nscphi/I"); my_tree->Branch("scphwh",scphwh,"scphwh[Nscphi]/I"); my_tree->Branch("scphstat",scphstat,"scphstat[Nscphi]/I"); my_tree->Branch("scphsect",scphsect,"scphsect[Nscphi]/I"); my_tree->Branch("scphbx",scphbx,"scphbx[Nscphi]/I"); my_tree->Branch("scphcod",scphcod,"scphcod[Nscphi]/I"); my_tree->Branch("scphphi",scphphi,"scphphi[Nscphi]/I"); my_tree->Branch("scphphib",scphphib,"scphphib[Nscphi]/I"); my_tree->Branch("scphposx",scphposx,"scphposx[Nscphi]/F"); my_tree->Branch("scphposy",scphposy,"scphposy[Nscphi]/F"); my_tree->Branch("scphposz",scphposz,"scphposz[Nscphi]/F"); my_tree->Branch("scphdirx",scphdirx,"scphdirx[Nscphi]/F"); my_tree->Branch("scphdiry",scphdiry,"scphdiry[Nscphi]/F"); my_tree->Branch("scphdirz",scphdirz,"scphdirz[Nscphi]/F"); // SC THETA block my_tree->Branch("Nsctheta",&nsctheta,"Nsctheta/I"); my_tree->Branch("scthwh",scthwh,"scthwh[Nsctheta]/I"); my_tree->Branch("scthstat",scthstat,"scthstat[Nsctheta]/I"); my_tree->Branch("scthsect",scthsect,"scthsect[Nsctheta]/I"); my_tree->Branch("scthbx",scthbx,"scthbx[Nsctheta]/I"); my_tree->Branch("scthcode",scthcode,"scthcode[Nsctheta][7]/I"); my_tree->Branch("scthpos",scthpos,"scthpos[Nsctheta][7]/I"); my_tree->Branch("scthqual",scthqual,"scthqual[Nsctheta][7]/I"); }
void DTTrigTest::beginRun | ( | const edm::Run & | iRun, |
const edm::EventSetup & | iEventSetup | ||
) | [virtual] |
Create DTTrig instance and TUs.
Reimplemented from edm::EDAnalyzer.
Definition at line 203 of file DTTrigTest.cc.
References gather_cfg::cout, DTTrig::createTUs(), my_debug, my_params, and my_trig.
void DTTrigTest::endJob | ( | void | ) | [virtual] |
Close Tree and write File.
Reimplemented from edm::EDAnalyzer.
Definition at line 76 of file DTTrigTest.cc.
References gather_cfg::cout, my_debug, my_rootfile, and my_tree.
{ if (my_debug) cout << "[DTTrigTest] Writing Tree and Closing File" << endl; my_tree->Write(); delete my_tree; my_rootfile->Close(); }
int DTTrigTest::bbx[100] [private] |
Definition at line 101 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::bcod[100] [private] |
Definition at line 102 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::bdirx[100] [private] |
Definition at line 108 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::bdiry[100] [private] |
Definition at line 109 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::bdirz[100] [private] |
Definition at line 110 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::bk[100] [private] |
Definition at line 103 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::bnum[100] [private] |
Definition at line 100 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::bposx[100] [private] |
Definition at line 105 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::bposy[100] [private] |
Definition at line 106 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::bposz[100] [private] |
Definition at line 107 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::bsect[100] [private] |
Definition at line 98 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::bsl[100] [private] |
Definition at line 99 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::bstat[100] [private] |
Definition at line 97 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::bwh[100] [private] |
Definition at line 96 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::bx[100] [private] |
Definition at line 104 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::chagen[10] [private] |
Definition at line 89 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::etagen[10] [private] |
Definition at line 87 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::eventn [private] |
Definition at line 78 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
bool DTTrigTest::my_debug [private] |
Definition at line 65 of file DTTrigTest.h.
Referenced by analyze(), beginRun(), DTTrigTest(), endJob(), and ~DTTrigTest().
edm::ParameterSet DTTrigTest::my_params [private] |
Definition at line 68 of file DTTrigTest.h.
Referenced by beginRun(), and DTTrigTest().
TFile* DTTrigTest::my_rootfile [private] |
Definition at line 73 of file DTTrigTest.h.
Referenced by DTTrigTest(), endJob(), and ~DTTrigTest().
TTree* DTTrigTest::my_tree [private] |
Definition at line 71 of file DTTrigTest.h.
Referenced by analyze(), beginJob(), DTTrigTest(), and endJob().
DTTrig* DTTrigTest::my_trig [private] |
Definition at line 62 of file DTTrigTest.h.
Referenced by analyze(), beginRun(), and ~DTTrigTest().
const double DTTrigTest::my_TtoTDC = 32./25. [static, private] |
Definition at line 59 of file DTTrigTest.h.
int DTTrigTest::nbti [private] |
Definition at line 95 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::ngen [private] |
Definition at line 82 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::nscphi [private] |
Definition at line 156 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::nsctheta [private] |
Definition at line 172 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::ntraco [private] |
Definition at line 113 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::ntsphi [private] |
Definition at line 130 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::ntstheta [private] |
Definition at line 146 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::phigen[10] [private] |
Definition at line 88 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::ptgen[10] [private] |
Definition at line 86 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::pxgen[10] [private] |
Definition at line 83 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::pygen[10] [private] |
Definition at line 84 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::pzgen[10] [private] |
Definition at line 85 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::runn [private] |
Definition at line 77 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::sbx[40] [private] |
Definition at line 134 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::scod[40] [private] |
Definition at line 135 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::scphbx[40] [private] |
Definition at line 160 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::scphcod[40] [private] |
Definition at line 161 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::scphdirx[100] [private] |
Definition at line 167 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::scphdiry[100] [private] |
Definition at line 168 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::scphdirz[100] [private] |
Definition at line 169 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::scphphi[40] [private] |
Definition at line 162 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::scphphib[40] [private] |
Definition at line 163 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::scphposx[100] [private] |
Definition at line 164 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::scphposy[100] [private] |
Definition at line 165 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::scphposz[100] [private] |
Definition at line 166 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::scphsect[40] [private] |
Definition at line 159 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::scphstat[40] [private] |
Definition at line 158 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::scphwh[40] [private] |
Definition at line 157 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::scthbx[40] [private] |
Definition at line 176 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::scthcode[40][7] [private] |
Definition at line 177 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::scthpos[40][7] [private] |
Definition at line 178 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::scthqual[40][7] [private] |
Definition at line 179 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::scthsect[40] [private] |
Definition at line 175 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::scthstat[40] [private] |
Definition at line 174 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::scthwh[40] [private] |
Definition at line 173 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::sdirx[100] [private] |
Definition at line 141 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::sdiry[100] [private] |
Definition at line 142 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::sdirz[100] [private] |
Definition at line 143 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::sphi[40] [private] |
Definition at line 136 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::sphib[40] [private] |
Definition at line 137 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::sposx[100] [private] |
Definition at line 138 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::sposy[100] [private] |
Definition at line 139 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::sposz[100] [private] |
Definition at line 140 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::ssect[40] [private] |
Definition at line 133 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::sstat[40] [private] |
Definition at line 132 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::swh[40] [private] |
Definition at line 131 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::tbx[80] [private] |
Definition at line 118 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::tcod[80] [private] |
Definition at line 119 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::tdirx[100] [private] |
Definition at line 125 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::tdiry[100] [private] |
Definition at line 126 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::tdirz[100] [private] |
Definition at line 127 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::thbx[40] [private] |
Definition at line 150 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::thcode[40][7] [private] |
Definition at line 151 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::thpos[40][7] [private] |
Definition at line 152 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::thqual[40][7] [private] |
Definition at line 153 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::thsect[40] [private] |
Definition at line 149 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::thstat[40] [private] |
Definition at line 148 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::thwh[40] [private] |
Definition at line 147 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::tk[80] [private] |
Definition at line 120 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::tnum[80] [private] |
Definition at line 117 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::tposx[100] [private] |
Definition at line 122 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::tposy[100] [private] |
Definition at line 123 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::tposz[100] [private] |
Definition at line 124 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::tsect[80] [private] |
Definition at line 116 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::tstat[80] [private] |
Definition at line 115 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::twh[80] [private] |
Definition at line 114 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
int DTTrigTest::tx[80] [private] |
Definition at line 121 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::vxgen[10] [private] |
Definition at line 90 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::vygen[10] [private] |
Definition at line 91 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::vzgen[10] [private] |
Definition at line 92 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().
float DTTrigTest::weight [private] |
Definition at line 79 of file DTTrigTest.h.
Referenced by analyze(), and beginJob().