Main Page   Namespace List   Class Hierarchy   Compound List   File List   Namespace Members   Compound Members   File Members  

Atlfast::JetMaker Class Reference

#include <JetMaker.h>

Collaboration diagram for Atlfast::JetMaker:

Collaboration graph
[legend]
List of all members.

Public Methods

 JetMaker (const std::string &name, ISvcLocator *pSvcLocator)
virtual ~JetMaker ()
StatusCode initialize ()
StatusCode execute ()
StatusCode finalize ()
Jet * create (Cluster *)
bool isAcceptable (Jet *)
void addMuons (Jet *, MuonCollection &)
MsgStream & dumpParameters (MsgStream &) const

Private Types

typedef std::vector< Cluster * > localClusterCollection
typedef localClusterCollection::iterator localClusterIterator
typedef std::vector< ReconstructedParticle * > MuonCollection
typedef MuonCollection::iterator MuonIterator

Private Methods

int getJetTags (MCparticleCollection &mcParticles, const PartitionCondition::BelowThresholdDeltaR &tagSelector) const
int getTauTags (MCparticleCollection &mcParticles, const PartitionCondition::BelowThresholdDeltaR &tagSelector, double jetPt) const
void tag (Jet *jet, MCparticleCollection &mcParticles_b, MCparticleCollection &mcParticles_c, MCparticleCollection &mcParticles_tau) const
HepLorentzVector addUnusedCells (MsgStream &) const

Private Attributes

double m_minPT
double m_maxEta
bool m_doSmearing
double m_barrelForwardEta
double m_rconeb
double m_rconef
double m_bPtMin
double m_bMaxDeltaR
double m_cPtMin
double m_cMaxDeltaR
double m_tauPtMin
double m_tauMaxDeltaR
double m_etaTagMax
double m_tauJetPtRatio
ISmearerm_smearer
ISmearerm_cellSmearer
TesIO * m_tesIO
std::string m_inputLocation
std::string m_outputLocation
std::string m_muonLocation
std::string m_mcTruthLocation
std::string m_unusedCellLocation
std::string m_missingMomentumLocation
KinematicHelper m_kinehelp

Member Typedef Documentation

typedef std::vector<Cluster*> Atlfast::JetMaker::localClusterCollection [private]
 

Definition at line 98 of file JetMaker.h.

Referenced by execute().

typedef localClusterCollection::iterator Atlfast::JetMaker::localClusterIterator [private]
 

Definition at line 99 of file JetMaker.h.

Referenced by execute().

typedef std::vector<ReconstructedParticle*> Atlfast::JetMaker::MuonCollection [private]
 

Definition at line 100 of file JetMaker.h.

Referenced by execute().

typedef MuonCollection::iterator Atlfast::JetMaker::MuonIterator [private]
 

Definition at line 101 of file JetMaker.h.

Referenced by addMuons().


Constructor & Destructor Documentation

Atlfast::JetMaker::JetMaker const std::string &    name,
ISvcLocator *    pSvcLocator
 

Definition at line 83 of file JetMaker.cxx.

00084     : Algorithm( name, pSvcLocator ){
00085     
00086     // Setting the parameter defaults.
00087     m_minPT            = 10;         // Min jet trans momentum
00088     m_maxEta           = 5;          //Max jet eta 
00089     m_doSmearing       = true;       //smearing is done
00090     m_rconeb           = 0.400;
00091     m_rconef           = 0.400;
00092     
00093     m_bPtMin           = 5.0;        //min pt for bquark to tag jet
00094     m_bMaxDeltaR       = 0.2;        //max dR for bquark to tag jet
00095     
00096     m_cPtMin           = 5.0;        //min pt for cquark to tag jet
00097     m_cMaxDeltaR       = 0.2;
00098     
00099     m_tauPtMin         = 10.0;
00100     m_tauMaxDeltaR     = 0.3;
00101     
00102     m_etaTagMax        = 2.5;
00103 
00104     m_tauJetPtRatio    = 0.9;
00105     
00106     
00107     // Default paths for entities in the TES
00108     m_inputLocation             = "/Event/AtlfastClusters";
00109     m_outputLocation            = "/Event/AtlfastJets";
00110     m_muonLocation              = "/Event/AtlfastNonIsolatedMuons";
00111     m_unusedCellLocation        = "/Event/AtlfastUnusedCells";
00112     m_missingMomentumLocation   = "/Event/AtlfastMissingMomentum";
00113     m_mcTruthLocation           = "/Event/McEventCollection";
00114     
00115     // This is how you declare the paramemters to Gaudi so that
00116     // they can be over-written via the job options file
00117     declareProperty( "MinimumPT",       m_minPT ) ;
00118     declareProperty( "MaximumEta",      m_maxEta ) ;
00119     declareProperty( "DoSmearing",      m_doSmearing ); 
00120     declareProperty( "RconeB",          m_rconeb );
00121     declareProperty( "RconeF",          m_rconef );
00122     declareProperty( "bPtMin",          m_bPtMin );
00123     declareProperty( "bMaxDeltaR",      m_bMaxDeltaR );
00124     declareProperty( "cPtMin",          m_cPtMin );
00125     declareProperty( "cMaxDeltaR",      m_cMaxDeltaR );
00126     declareProperty( "tauPtMin",        m_tauPtMin );
00127     declareProperty( "tauMaxDeltaR",    m_tauMaxDeltaR );
00128     declareProperty( "etaTagMax",       m_etaTagMax );
00129     declareProperty( "tauJetPtRatio",   m_tauJetPtRatio);
00130     
00131     declareProperty( "InputLocation",           m_inputLocation ) ;
00132     declareProperty( "OutputLocation",          m_outputLocation ) ;
00133     declareProperty( "UnusedCellLocation",      m_unusedCellLocation ) ;
00134     declareProperty( "MuonLocation",            m_muonLocation );
00135     declareProperty( "MissingMomentumLocation", m_missingMomentumLocation );
00136     declareProperty( "McTruthLocation",         m_mcTruthLocation );
00137   }

Atlfast::JetMaker::~JetMaker   [virtual]
 

Definition at line 139 of file JetMaker.cxx.

References m_cellSmearer, m_smearer, and m_tesIO.

00139                       {
00140     
00141     MsgStream log( messageService(), name() ) ;
00142     log << MSG::INFO << "Destructor Called" << endreq;
00143     
00144     if (m_smearer) {
00145       log << MSG::INFO << "Deleting jet Smearer" << endreq;
00146       delete m_smearer;
00147     }
00148     if (m_cellSmearer) {
00149       log << MSG::INFO << "Deleting Cell Smearer" << endreq;
00150       delete m_cellSmearer;
00151     }
00152     if (m_tesIO) {
00153       delete m_tesIO;
00154     }
00155   } 

Member Function Documentation

int Atlfast::JetMaker::getJetTags MCparticleCollection &    mcParticles,
const PartitionCondition::BelowThresholdDeltaR &    tagSelector
const [private]
 

Definition at line 546 of file JetMaker.cxx.

00547                                                                                           {
00548     
00549     MsgStream log( messageService(), name() ) ;
00550 
00551     MCparticleCollectionCIter src ;
00552         for( src = mcParticles.begin(); src != mcParticles.end(); ++src ){
00553           if(tagSelector(SimpleKinematic(*src))) return (*src)->pdg_id();
00554         }
00555         return 0;
00556   }

int Atlfast::JetMaker::getTauTags MCparticleCollection &    mcParticles,
const PartitionCondition::BelowThresholdDeltaR &    tagSelector,
double    jetPt
const [private]
 

Definition at line 560 of file JetMaker.cxx.

00564                                                    {
00565     
00566     MsgStream log( messageService(), name() ) ;
00567     MCparticleCollectionCIter src ;
00568     for( src = mcParticles.begin() ; src != mcParticles.end() ; ++src )
00569       {
00570         double tauPt = (*src)->momentum().perp();
00571         double neuPt=0.;
00572         // Find Tau children
00573         HepMC::GenVertex*  endVertex = (*src)->end_vertex();
00574         if(!endVertex) return 0;
00575 
00576         HepMC::GenVertex::particle_iterator firstChild = 
00577           endVertex->particles_begin(HepMC::children);
00578         
00579         HepMC::GenVertex::particle_iterator endChild = 
00580           endVertex->particles_end(HepMC::children);
00581         
00582         HepMC::GenVertex::particle_iterator thisChild = firstChild; 
00583         
00584         // Find Pt of Taus tau-neutrino
00585         for(; thisChild!=endChild; ++thisChild){
00586           if(  abs((*thisChild)->pdg_id()) == 16){
00587             neuPt = (*thisChild)->momentum().perp();
00588           }
00589         }
00590           
00591           if(tagSelector(SimpleKinematic(*src))      &&
00592              ((tauPt-neuPt)/jetPt) >= m_tauJetPtRatio) return (*src)->pdg_id();
00593         
00594         }
00595       return 0;
00596   }   

void Atlfast::JetMaker::tag Jet *    jet,
MCparticleCollection &    mcParticles_b,
MCparticleCollection &    mcParticles_c,
MCparticleCollection &    mcParticles_tau
const [private]
 

Definition at line 512 of file JetMaker.cxx.

Referenced by execute().

00515                                                                  {
00516  
00517     //See if jet can be b-tagged
00518     //Construct an object to select IKinematics within DeltaR of the current jet
00519     PartitionCondition::BelowThresholdDeltaR bTagSelector(jet, m_bMaxDeltaR); 
00520     int btag = getJetTags(mcParticles_b,bTagSelector);
00521     if(btag != 0 )
00522       {
00523         jet->setBTag();
00524         return;
00525       }
00526     //else see if jet can be c-tagged
00527     //Construct an object to select IKinematics within DeltaR of the current jet
00528     PartitionCondition::BelowThresholdDeltaR cTagSelector(jet, m_cMaxDeltaR); 
00529     int ctag = getJetTags( mcParticles_c, cTagSelector);
00530     if(ctag != 0 )
00531       { 
00532         jet->setCTag();
00533         return;
00534       }
00535 
00536     //finally see if jet can be tau-tagged
00537     //Construct an object to select IKinematics within DeltaR of the current jet
00538     PartitionCondition::BelowThresholdDeltaR tauTagSelector(jet, m_tauMaxDeltaR); 
00539     int tauTag=getTauTags( mcParticles_tau, tauTagSelector, jet->pT());
00540     if(tauTag !=0 ) jet->setTauTag(tauTag);
00541     return; 
00542   }

HepLorentzVector Atlfast::JetMaker::addUnusedCells MsgStream &    const [private]
 

Definition at line 598 of file JetMaker.cxx.

References m_cellSmearer, m_tesIO, and m_unusedCellLocation.

Referenced by execute().

00598                                                                {
00599 
00600     HepLorentzVector temp(0., 0., 0., 0.);
00601     
00602 
00603     DataHandle<ITwoCptCellCollection> unusedCells;
00604     if( ! m_tesIO->getDH( unusedCells, m_unusedCellLocation) ){
00605       log << MSG::INFO << "No unused cells  in TES " << endreq;
00606       return temp;
00607     }  
00608     log << MSG::DEBUG << "No unused cells  in TES " << endreq;
00609     
00610     return std::accumulate(unusedCells->begin(), 
00611                            unusedCells->end(), 
00612                            temp, 
00613                            SmearCell(m_cellSmearer)
00614                            );
00615 
00616   }

StatusCode Atlfast::JetMaker::initialize  
 

Definition at line 161 of file JetMaker.cxx.

References Atlfast::GlobalEventData::barrelForwardEta(), Atlfast::GlobalEventData::lumi(), m_barrelForwardEta, m_bMaxDeltaR, m_bPtMin, m_cellSmearer, m_cMaxDeltaR, m_cPtMin, m_doSmearing, m_etaTagMax, m_inputLocation, m_maxEta, m_mcTruthLocation, m_minPT, m_missingMomentumLocation, m_muonLocation, m_outputLocation, m_rconeb, m_rconef, m_smearer, m_tauJetPtRatio, m_tauMaxDeltaR, m_tauPtMin, m_tesIO, and Atlfast::GlobalEventData::randSeed().

00161                                  {
00162     MsgStream log( messageService(), name() ) ;
00163     
00164     log << MSG::DEBUG << "instantiating a JetSmearer" << endreq;
00165     
00166     //    m_tesIO = new TesIO(eventDataService());  
00167     m_tesIO = new TesIO();  
00168     
00169     //get the Global Event Data using singleton pattern
00170     GlobalEventData* ged = GlobalEventData::Instance();
00171     int lumi                  = ged->lumi();
00172     int randSeed              = ged->randSeed() ;
00173     m_barrelForwardEta = ged->barrelForwardEta();
00174     
00175     if(m_doSmearing){
00176       m_smearer                 = new JetSmearer(randSeed, 
00177                                                  lumi,
00178                                                  m_rconeb, 
00179                                                  m_rconef, 
00180                                                  m_barrelForwardEta);    
00181       m_cellSmearer             = new JetSmearer(randSeed, 
00182                                                  lumi,
00183                                                  0., 
00184                                                  0., 
00185                                                  m_barrelForwardEta);    
00186     } else {
00187       m_smearer = NULL;
00188       m_cellSmearer = NULL;
00189     }
00190 
00191     HeaderPrinter hp("Atlfast JetMaker:", log);
00192     
00193     hp.add("TES Locations:              ");
00194     hp.add(" Muons from                 ", m_muonLocation); 
00195     hp.add(" Clusters from              ", m_inputLocation); 
00196     hp.add(" GeneratorParticles from    ", m_mcTruthLocation);
00197     hp.add(" Jets to                    ", m_outputLocation);
00198     hp.add(" unused cell,cluster sum to ", m_missingMomentumLocation);
00199     hp.add("Smearing on                 ", m_doSmearing);
00200     hp.add("Cluster min Pt              ", m_minPT);
00201     hp.add("Cluster max Eta             ", m_maxEta);           
00202     hp.add("Luminosity                  ", lumi);
00203     hp.add("Initial random seed         ", randSeed);      
00204     hp.add("End cap Cone size           ", m_rconef);
00205     hp.add("Barrel  Cone size           ", m_rconeb);
00206     hp.add("Eta of start of end cap     ", m_barrelForwardEta);
00207     hp.add("Parameters for tags         ", "       ");
00208     hp.add(" Max eta                    ", m_etaTagMax);
00209     hp.add(" b-quark min pt             ", m_bPtMin);
00210     hp.add(" b-quark max DeltaR         ", m_bMaxDeltaR);
00211     hp.add(" c-quark min pt             ", m_cPtMin);
00212     hp.add(" c-quark max DeltaR         ", m_cMaxDeltaR);
00213     hp.add(" tau lep min pt             ", m_tauPtMin);
00214     hp.add(" tau lep max DeltaR         ", m_tauMaxDeltaR);
00215     hp.add(" tau/Jet Pt Ratio           ", m_tauJetPtRatio);
00216     hp.print();
00217     
00218 
00219 
00220     return StatusCode::SUCCESS ;
00221   }

StatusCode Atlfast::JetMaker::execute  
 

Definition at line 252 of file JetMaker.cxx.

References addMuons(), addUnusedCells(), create(), isAcceptable(), localClusterCollection, localClusterIterator, m_bPtMin, m_cPtMin, m_etaTagMax, m_inputLocation, m_missingMomentumLocation, m_muonLocation, m_outputLocation, m_tauPtMin, m_tesIO, MuonCollection, and tag().

00252                                {
00253     
00254     MsgStream log( messageService(), name() ) ;    
00255 
00256     std::string mess;
00257     // make a message logging stream
00258     
00259 
00260       
00261 
00262     //.............................
00263     // Make some locals stores which be used to keep pointers to
00264     // all of the entities from the event store which are needed by
00265     // this algorithm. These are all local and are typedefed in this class
00266     
00267     localClusterCollection      my_Clusters ;
00268     MuonCollection         my_Muons ;
00269     //................................
00270     
00271     //    if( ! m_tesIO->copy( my_Clusters, m_inputLocation ) ) {
00272     if( ! m_tesIO->copy<ClusterCollection> ( my_Clusters, m_inputLocation ) ) {
00273       log << MSG::INFO << "No Clusters in TES " << endreq;
00274     } else {
00275       log << MSG::DEBUG << my_Clusters.size()<<" Clusters from TES " << endreq;
00276     }
00277     // get non-isolated muons for jet construction
00278     
00279     if( ! m_tesIO->copy<ReconstructedParticleCollection>
00280         ( my_Muons, m_muonLocation ) ) {
00281       log << MSG::DEBUG << "No muons in TES " << endreq;
00282     }  
00283 
00284     // ......................
00285     // Make a container object which will store pointers to all isolated 
00286     // Jetss which are successfully created.  Since it is going to go 
00287     // into the event store, it has to be a special Gaudi type of container. 
00288     // This is typedefined in the entity include file
00289     // As it is going into the TES it needs to be created on the heap with new.
00290     
00291     JetCollection* myJets = new JetCollection ;
00292     
00293         //.........................................................
00294     // From each Cluster create a reconstructed Jet. 
00295     // If all requirements are satisfied add to the collection
00296     
00297     Jet* candidate ;
00298     localClusterIterator src ;
00299     HepLorentzVector missingMomentum(0.,0.,0.,0.);
00300     const int BQUARK(ParticleCodes::BQUARK);
00301     const int CQUARK(ParticleCodes::CQUARK);
00302     //Construct an object to select b, c and tau quarks with kine conditions
00303     const HepMC_helper::SelectJetTag bSelector(BQUARK, m_bPtMin, m_etaTagMax);
00304     const HepMC_helper::SelectJetTag cSelector(CQUARK, m_cPtMin, m_etaTagMax);
00305     const HepMC_helper::SelectTauTag tauSelector(m_tauPtMin, m_etaTagMax);
00306 
00307     MCparticleCollection  mcParticles_b ;
00308     MCparticleCollection  mcParticles_c ;
00309     MCparticleCollection  mcParticles_tau ;
00310     //get quarks used in tagging from TES using tesio
00311     m_tesIO->getMC( mcParticles_b, &bSelector ) ;
00312     m_tesIO->getMC( mcParticles_c, &cSelector ) ;
00313     m_tesIO->getMC( mcParticles_tau, &tauSelector ) ;
00314  
00315     for( src = my_Clusters.begin() ; src != my_Clusters.end() ; ++src ) 
00316       { 
00317         IAssociationManager* iAssocMan = *src;
00318         ReconstructedParticle rp;
00319         if (iAssocMan->associations(rp).size() == 0)
00320           {
00321             // Make a new jet
00322             candidate = this->create( *src );
00323 
00324             //Make a temp HepLorentzVector from candidates 4-mom
00325             //before muon is added in
00326             HepLorentzVector temp4Vec = candidate->momentum();
00327 
00328             // add isolated muons
00329             this->addMuons(candidate, my_Muons);
00330 
00331             
00332             // add to missing momentum vector if not used in a jet.
00333             // remember, these are smeared momenta if smearing is enabled
00334             if( this->isAcceptable( candidate ) ) {    
00335 
00336             // if the jet is acceptable flavor tag the jet
00337               this->tag(candidate,mcParticles_b,mcParticles_c,mcParticles_tau);
00338 
00339               log << MSG::DEBUG << "Pushing back a Jet" << endreq;
00340               myJets->push_back( candidate ) ; 
00341               log<<MSG::DEBUG<<"jet "<<*candidate<<endreq;
00342             }else {
00343               log << MSG::DEBUG << "Jet fails cuts" << endreq;
00344               //add in temp4Vec to avoid muon double counting
00345               missingMomentum += temp4Vec;
00346               delete candidate ;
00347             }
00348           }
00349       }
00350     
00351     log<<MSG::DEBUG<<"Missing momentum: clusters "<<missingMomentum<<endreq;
00352     missingMomentum += addUnusedCells(log);
00353     log<<MSG::DEBUG<<"Missing momentum: cells+clstrs "
00354        <<missingMomentum<<endreq;
00355     //......................................
00356     // Register any Jets which have been successfilly created in the 
00357     // transient event store. 
00358     TesIoStat stat;
00359     log<<MSG::DEBUG<<"Storing "<< myJets->size() <<" Jets"<<endreq;
00360     stat = m_tesIO->store( myJets, m_outputLocation ) ;
00361     if(!stat){
00362       log<<MSG::ERROR<<"Could not store jets"<<endreq;
00363       return stat;
00364     }
00365 
00366 
00367     MissingMomentum* mm = new MissingMomentum(missingMomentum);
00368     stat = m_tesIO->store(mm,m_missingMomentumLocation ) ;
00369     std::string message = stat? "Stored missing mom":"Failed to missing mom";
00370     log<<MSG::DEBUG<<message<<endreq;
00371   
00372 
00373     return StatusCode::SUCCESS;
00374   }

StatusCode Atlfast::JetMaker::finalize  
 

Definition at line 227 of file JetMaker.cxx.

00227                                {
00228 
00229     MsgStream log( messageService(), name() ) ;
00230     log << MSG::INFO << "finalizing" << endreq;    
00231     return StatusCode::SUCCESS ;
00232   }

Jet * Atlfast::JetMaker::create Cluster *   
 

Definition at line 389 of file JetMaker.cxx.

References m_smearer, and Atlfast::ISmearer::smear().

Referenced by execute().

00390   {
00391     Jet* jet;
00392     MsgStream log( messageService(), name() ) ;
00393     HepLorentzVector vec = src->momentum();
00394     
00395     if (m_doSmearing) 
00396       {
00397         // Ask our Smearer make a smeared 4-vector
00398         vec = m_smearer->smear(vec);      
00399       }
00400     
00401     jet = new Jet(vec, *src);
00402     log << MSG::DEBUG << "Unsmeared Cluster: " << *src << endreq;
00403     log << MSG::DEBUG << "Smeared Jet      : " << jet << endreq ;
00404     
00405     // return Jet
00406     return jet;
00407     
00408     
00409   }

bool Atlfast::JetMaker::isAcceptable Jet *   
 

Definition at line 484 of file JetMaker.cxx.

References m_maxEta, and m_minPT.

Referenced by execute().

00485   {
00486     
00487     MsgStream log( messageService(), name() ) ;
00488     
00489     // Apply kinematic criteria to the candidate jet
00490     
00491     if( candidate->pT() < m_minPT ) {        
00492       log << MSG::DEBUG << "Candidate failed pt cut: pt " << candidate->pT() 
00493           << " min was " << m_minPT << endreq;
00494       return false ;
00495     }
00496     
00497     if( abs(candidate->eta()) > m_maxEta ) {
00498       log << MSG::DEBUG << "Candidate failed max eta cut: eta " 
00499           << candidate->eta() 
00500           << " max was " << m_maxEta << endreq;
00501       return false ;
00502     }
00503     
00504     return true ;
00505   }

void Atlfast::JetMaker::addMuons Jet *   ,
MuonCollection  
 

Definition at line 416 of file JetMaker.cxx.

References m_barrelForwardEta, m_kinehelp, m_rconeb, m_rconef, and MuonIterator.

Referenced by execute().

00417   {
00418     
00419     if (muons.size() == 0)  return;
00420     
00421     MsgStream log( messageService(), name() ) ;
00422     
00423     MuonIterator first = muons.begin() ;
00424     MuonIterator last  = muons.end() ;
00425     
00426     // Sort by distance from the jet
00427     sort( first, last, SortAttribute::DeltaR(jet) );
00428     
00429     // check that closest muon is within cone  
00430     // if it is, add muon to jet
00431     
00432     float muonJetRdist = m_kinehelp.deltaR( jet, *first ) ;
00433     
00434     log << MSG::DEBUG << "distance in R of muon to Jet " << muonJetRdist << endreq;
00435     
00436     if (
00437         (
00438          abs((*first)->eta()) <= m_barrelForwardEta && 
00439          muonJetRdist < m_rconeb
00440          ) 
00441         ||
00442         (
00443          abs((*first)->eta()) >  m_barrelForwardEta && 
00444          muonJetRdist < m_rconef
00445          )
00446         ) 
00447       {
00448         
00449         log << MSG::DEBUG 
00450             << "Adding muon to jet with momentum "
00451             << (*first)->momentum() 
00452             << endreq;
00453         
00454         log << MSG::DEBUG << "Old mom: " << *jet << endreq;
00455         
00456         jet->setMomentum(jet->momentum() + (*first)->momentum()); 
00457         
00458         log << MSG::DEBUG << "New mom: " << *jet << endreq ;
00459         
00460         // create association
00461         IAssociationManager* ia = jet;
00462         ia->associate(*first);
00463         
00464         // remove muon from list
00465         muons.erase(first);
00466         
00467         
00468       }
00469     else
00470       {        
00471         log << MSG::DEBUG << "Muon not in R cone" << endreq;
00472       }
00473     
00474   }

MsgStream& Atlfast::JetMaker::dumpParameters MsgStream &    const
 


Member Data Documentation

double Atlfast::JetMaker::m_minPT [private]
 

Definition at line 109 of file JetMaker.h.

Referenced by initialize(), and isAcceptable().

double Atlfast::JetMaker::m_maxEta [private]
 

Definition at line 110 of file JetMaker.h.

Referenced by initialize(), and isAcceptable().

bool Atlfast::JetMaker::m_doSmearing [private]
 

Definition at line 111 of file JetMaker.h.

Referenced by initialize().

double Atlfast::JetMaker::m_barrelForwardEta [private]
 

Definition at line 113 of file JetMaker.h.

Referenced by addMuons(), and initialize().

double Atlfast::JetMaker::m_rconeb [private]
 

Definition at line 114 of file JetMaker.h.

Referenced by addMuons(), and initialize().

double Atlfast::JetMaker::m_rconef [private]
 

Definition at line 115 of file JetMaker.h.

Referenced by addMuons(), and initialize().

double Atlfast::JetMaker::m_bPtMin [private]
 

Definition at line 117 of file JetMaker.h.

Referenced by execute(), and initialize().

double Atlfast::JetMaker::m_bMaxDeltaR [private]
 

Definition at line 118 of file JetMaker.h.

Referenced by initialize().

double Atlfast::JetMaker::m_cPtMin [private]
 

Definition at line 119 of file JetMaker.h.

Referenced by execute(), and initialize().

double Atlfast::JetMaker::m_cMaxDeltaR [private]
 

Definition at line 120 of file JetMaker.h.

Referenced by initialize().

double Atlfast::JetMaker::m_tauPtMin [private]
 

Definition at line 121 of file JetMaker.h.

Referenced by execute(), and initialize().

double Atlfast::JetMaker::m_tauMaxDeltaR [private]
 

Definition at line 122 of file JetMaker.h.

Referenced by initialize().

double Atlfast::JetMaker::m_etaTagMax [private]
 

Definition at line 123 of file JetMaker.h.

Referenced by execute(), and initialize().

double Atlfast::JetMaker::m_tauJetPtRatio [private]
 

Definition at line 124 of file JetMaker.h.

Referenced by initialize().

ISmearer* Atlfast::JetMaker::m_smearer [private]
 

Definition at line 132 of file JetMaker.h.

Referenced by create(), initialize(), and ~JetMaker().

ISmearer* Atlfast::JetMaker::m_cellSmearer [private]
 

Definition at line 133 of file JetMaker.h.

Referenced by addUnusedCells(), initialize(), and ~JetMaker().

TesIO* Atlfast::JetMaker::m_tesIO [private]
 

Definition at line 136 of file JetMaker.h.

Referenced by addUnusedCells(), execute(), initialize(), and ~JetMaker().

std::string Atlfast::JetMaker::m_inputLocation [private]
 

Definition at line 144 of file JetMaker.h.

Referenced by execute(), and initialize().

std::string Atlfast::JetMaker::m_outputLocation [private]
 

Definition at line 145 of file JetMaker.h.

Referenced by execute(), and initialize().

std::string Atlfast::JetMaker::m_muonLocation [private]
 

Definition at line 146 of file JetMaker.h.

Referenced by execute(), and initialize().

std::string Atlfast::JetMaker::m_mcTruthLocation [private]
 

Definition at line 147 of file JetMaker.h.

Referenced by initialize().

std::string Atlfast::JetMaker::m_unusedCellLocation [private]
 

Definition at line 148 of file JetMaker.h.

Referenced by addUnusedCells().

std::string Atlfast::JetMaker::m_missingMomentumLocation [private]
 

Definition at line 149 of file JetMaker.h.

Referenced by execute(), and initialize().

KinematicHelper Atlfast::JetMaker::m_kinehelp [private]
 

Definition at line 155 of file JetMaker.h.

Referenced by addMuons().


The documentation for this class was generated from the following files:
Generated on Tue Mar 18 11:18:56 2003 for AtlfastAlgs by doxygen1.3-rc1