#include <RecoLocalCalo/Castor/src/CastorTowerProducer.cc>
Public Member Functions | |
CastorTowerProducer (const edm::ParameterSet &) | |
~CastorTowerProducer () | |
Private Types | |
typedef edm::SortedCollection < CastorRecHit > | CastorRecHitCollection |
typedef edm::RefVector < CastorRecHitCollection > | CastorRecHitRefVector |
typedef std::vector < reco::CastorTower > | CastorTowerCollection |
typedef ROOT::Math::RhoZPhiPoint | CellPoint |
typedef math::XYZPointD | Point |
typedef ROOT::Math::RhoEtaPhiPoint | TowerPoint |
Private Member Functions | |
virtual void | beginJob () |
virtual void | ComputeTowerVariable (const edm::RefVector< edm::SortedCollection< CastorRecHit > > &usedRecHits, double &Ehot, double &depth) |
virtual void | endJob () |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
Private Attributes | |
std::string | input_ |
double | maxtime_ |
double | mintime_ |
double | towercut_ |
Description: CastorTower Reconstruction Producer. Produce CastorTowers from CastorCells. Implementation:
Definition at line 45 of file CastorTowerProducer.cc.
typedef edm::SortedCollection<CastorRecHit> CastorTowerProducer::CastorRecHitCollection [private] |
Definition at line 60 of file CastorTowerProducer.cc.
typedef edm::RefVector<CastorRecHitCollection> CastorTowerProducer::CastorRecHitRefVector [private] |
Definition at line 62 of file CastorTowerProducer.cc.
typedef std::vector<reco::CastorTower> CastorTowerProducer::CastorTowerCollection [private] |
Definition at line 61 of file CastorTowerProducer.cc.
typedef ROOT::Math::RhoZPhiPoint CastorTowerProducer::CellPoint [private] |
Definition at line 59 of file CastorTowerProducer.cc.
typedef math::XYZPointD CastorTowerProducer::Point [private] |
Definition at line 57 of file CastorTowerProducer.cc.
typedef ROOT::Math::RhoEtaPhiPoint CastorTowerProducer::TowerPoint [private] |
Definition at line 58 of file CastorTowerProducer.cc.
CastorTowerProducer::CastorTowerProducer | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 83 of file CastorTowerProducer.cc.
: input_(iConfig.getUntrackedParameter<std::string>("inputprocess","castorreco")), towercut_(iConfig.getUntrackedParameter<double>("towercut",1.)), mintime_(iConfig.getUntrackedParameter<double>("mintime",-999)), maxtime_(iConfig.getUntrackedParameter<double>("maxtime",999)) { //register your products produces<CastorTowerCollection>(); //now do what ever other initialization is needed }
CastorTowerProducer::~CastorTowerProducer | ( | ) |
Definition at line 95 of file CastorTowerProducer.cc.
{ // do anything here that needs to be done at desctruction time // (e.g. close files, deallocate resources etc.) }
void CastorTowerProducer::beginJob | ( | void | ) | [private, virtual] |
Reimplemented from edm::EDProducer.
Definition at line 247 of file CastorTowerProducer.cc.
References LogDebug.
{ LogDebug("CastorTowerProducer") <<"Starting CastorTowerProducer"; }
void CastorTowerProducer::ComputeTowerVariable | ( | const edm::RefVector< edm::SortedCollection< CastorRecHit > > & | usedRecHits, |
double & | Ehot, | ||
double & | depth | ||
) | [private, virtual] |
Definition at line 258 of file CastorTowerProducer.cc.
References edm::Ref< C, T, F >::id(), module::module(), and dt_dqm_sourceclient_common_cff::reco.
Referenced by produce().
{ using namespace reco; double Etot = 0; // loop over the cells used in the tower k for (CastorRecHitRefVector::iterator it = usedRecHits.begin(); it != usedRecHits.end(); it++) { edm::Ref<CastorRecHitCollection> rechit_p = *it; double Erechit = rechit_p->energy(); HcalCastorDetId id = rechit_p->id(); int module = id.module(); double zrechit = 0; if (module < 3) zrechit = -14390 - 24.75 - 49.5*(module-1); if (module > 2) zrechit = -14390 - 99 - 49.5 - 99*(module-3); if(Erechit > Ehot) Ehot = Erechit; depth+=Erechit*zrechit; Etot+=Erechit; } depth/=Etot; Ehot/=Etot; }
void CastorTowerProducer::endJob | ( | void | ) | [private, virtual] |
Reimplemented from edm::EDProducer.
Definition at line 253 of file CastorTowerProducer.cc.
References LogDebug.
{ LogDebug("CastorTowerProducer") <<"Ending CastorTowerProducer"; }
void CastorTowerProducer::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Implements edm::EDProducer.
Definition at line 107 of file CastorTowerProducer.cc.
References ComputeTowerVariable(), edm::Event::getByLabel(), i, edm::Ref< C, T, F >::id(), input_, j, gen::k, LogDebug, maxtime_, mintime_, module::module(), edm::Event::put(), dt_dqm_sourceclient_common_cff::reco, and towercut_.
{ using namespace edm; using namespace reco; using namespace TMath; // Produce CastorTowers from CastorCells edm::Handle<CastorRecHitCollection> InputRecHits; iEvent.getByLabel(input_,InputRecHits); std::auto_ptr<CastorTowerCollection> OutputTowers (new CastorTowerCollection); // get and check input size int nRecHits = InputRecHits->size(); LogDebug("CastorTowerProducer") <<"2. entering CastorTowerProducer"<<std::endl; if (nRecHits==0) LogDebug("CastorTowerProducer") <<"Warning: You are trying to run the Tower algorithm with 0 input rechits."; // declare castor array // (0,x): Energies - (1,x): emEnergies - (2,x): hadEnergies - (3,x): phi position double poscastortowerarray[4][16]; double negcastortowerarray[4][16]; CastorRecHitRefVector poscastorusedrechits[16]; CastorRecHitRefVector negcastorusedrechits[16]; // set phi values and everything else to zero for (int j = 0; j < 16; j++) { poscastortowerarray[3][j] = -2.94524 + j*0.3927; poscastortowerarray[0][j] = 0.; poscastortowerarray[1][j] = 0.; poscastortowerarray[2][j] = 0.; negcastortowerarray[3][j] = -2.94524 + j*0.3927; negcastortowerarray[0][j] = 0.; negcastortowerarray[1][j] = 0.; negcastortowerarray[2][j] = 0.; } // loop over rechits to build castortowerarray[4][16] and castorusedrechits[16] for (unsigned int i = 0; i < InputRecHits->size(); i++) { edm::Ref<CastorRecHitCollection> rechit_p = edm::Ref<CastorRecHitCollection>(InputRecHits, i); double Erechit = rechit_p->energy(); HcalCastorDetId id = rechit_p->id(); int module = id.module(); int sector = id.sector(); double zrechit = 0; if (module < 3) zrechit = -14390 - 24.75 - 49.5*(module-1); if (module > 2) zrechit = -14390 - 99 - 49.5 - 99*(module-3); double phirechit = -100; if (sector < 9) phirechit = 0.19635 + (sector-1)*0.3927; if (sector > 8) phirechit = -2.94524 + (sector - 9)*0.3927; // add time conditions for the rechit if (rechit_p->time() > mintime_ && rechit_p->time() < maxtime_) { // loop over the 16 towers possibilities for ( int j=0;j<16;j++) { // phi matching condition if (TMath::Abs(phirechit - poscastortowerarray[3][j]) < 0.0001) { // condition over rechit z value if (zrechit > 0.) { poscastortowerarray[0][j]+=Erechit; if (module < 3) {poscastortowerarray[1][j]+=Erechit;} else {poscastortowerarray[2][j]+=Erechit;} poscastorusedrechits[j].push_back(rechit_p); } else { negcastortowerarray[0][j]+=Erechit; if (module < 3) {negcastortowerarray[1][j]+=Erechit;} else {negcastortowerarray[2][j]+=Erechit;} negcastorusedrechits[j].push_back(rechit_p); } // end condition over rechit z value } // end phi matching condition } // end loop over the 16 towers possibilities } // end time conditions } // end loop over rechits to build castortowerarray[4][16] and castorusedrechits[16] // make towers of the arrays double fem, Ehot, depth; double rhoTower = 88.5; // loop over the 16 towers possibilities for (int k=0;k<16;k++) { fem = 0; Ehot = 0; depth = 0; // select the positive towers with E > Ecut if (poscastortowerarray[0][k] > towercut_) { fem = poscastortowerarray[1][k]/poscastortowerarray[0][k]; CastorRecHitRefVector usedRecHits = poscastorusedrechits[k]; ComputeTowerVariable(usedRecHits,Ehot,depth); LogDebug("CastorTowerProducer") <<"tower "<<k+1<<": fem = "<<fem<<" ,depth = "<<depth<<" ,Ehot = "<<Ehot<<std::endl; TowerPoint temptowerposition(rhoTower,5.9,poscastortowerarray[3][k]); Point towerposition(temptowerposition); CastorTower newtower(poscastortowerarray[0][k],towerposition,poscastortowerarray[1][k],poscastortowerarray[2][k],fem,depth,Ehot, poscastorusedrechits[k]); OutputTowers->push_back(newtower); } // end select the positive towers with E > Ecut // select the negative towers with E > Ecut if (negcastortowerarray[0][k] > towercut_) { fem = negcastortowerarray[1][k]/negcastortowerarray[0][k]; CastorRecHitRefVector usedRecHits = negcastorusedrechits[k]; ComputeTowerVariable(usedRecHits,Ehot,depth); LogDebug("CastorTowerProducer") <<"tower "<<k+1 << " energy = " << negcastortowerarray[0][k] << "EM = " << negcastortowerarray[1][k] << "HAD = " << negcastortowerarray[2][k] << "phi = " << negcastortowerarray[3][k] << ": fem = "<<fem<<" ,depth = "<<depth<<" ,Ehot = "<<Ehot<<std::endl; TowerPoint temptowerposition(rhoTower,-5.9,negcastortowerarray[3][k]); Point towerposition(temptowerposition); CastorTower newtower(negcastortowerarray[0][k],towerposition,negcastortowerarray[1][k],negcastortowerarray[2][k],fem,depth,Ehot, negcastorusedrechits[k]); OutputTowers->push_back(newtower); } // end select the negative towers with E > Ecut } // end loop over the 16 towers possibilities iEvent.put(OutputTowers); }
std::string CastorTowerProducer::input_ [private] |
Definition at line 63 of file CastorTowerProducer.cc.
Referenced by produce().
double CastorTowerProducer::maxtime_ [private] |
Definition at line 66 of file CastorTowerProducer.cc.
Referenced by produce().
double CastorTowerProducer::mintime_ [private] |
Definition at line 65 of file CastorTowerProducer.cc.
Referenced by produce().
double CastorTowerProducer::towercut_ [private] |
Definition at line 64 of file CastorTowerProducer.cc.
Referenced by produce().