CMS 3D CMS Logo

Public Member Functions | Private Types | Private Member Functions | Private Attributes

CastorClusterProducer Class Reference

#include <RecoLocalCalo/Castor/src/CastorClusterProducer.cc>

Inheritance diagram for CastorClusterProducer:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 CastorClusterProducer (const edm::ParameterSet &)
 ~CastorClusterProducer ()

Private Types

typedef std::vector
< reco::CastorCluster
CastorClusterCollection
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 endJob ()
double phiangle (double testphi)
virtual void produce (edm::Event &, const edm::EventSetup &)

Private Attributes

std::string basicjets_
bool clusteralgo_
std::string input_

Detailed Description

Description: CastorCluster Reconstruction Producer. Produces Clusters from Towers Implementation:

Definition at line 50 of file CastorClusterProducer.cc.


Member Typedef Documentation

Definition at line 66 of file CastorClusterProducer.cc.

Definition at line 65 of file CastorClusterProducer.cc.

typedef ROOT::Math::RhoZPhiPoint CastorClusterProducer::CellPoint [private]

Definition at line 64 of file CastorClusterProducer.cc.

Definition at line 62 of file CastorClusterProducer.cc.

typedef ROOT::Math::RhoEtaPhiPoint CastorClusterProducer::TowerPoint [private]

Definition at line 63 of file CastorClusterProducer.cc.


Constructor & Destructor Documentation

CastorClusterProducer::CastorClusterProducer ( const edm::ParameterSet iConfig) [explicit]

Definition at line 85 of file CastorClusterProducer.cc.

                                                                           :
  input_(iConfig.getUntrackedParameter<std::string>("inputtowers","")),
  basicjets_(iConfig.getUntrackedParameter<std::string>("basicjets","")),
  clusteralgo_(iConfig.getUntrackedParameter<bool>("ClusterAlgo",false))
{
  // register your products
  produces<CastorClusterCollection>();
  // now do what ever other initialization is needed
}
CastorClusterProducer::~CastorClusterProducer ( )

Definition at line 96 of file CastorClusterProducer.cc.

{
  // do anything here that needs to be done at desctruction time
  // (e.g. close files, deallocate resources etc.)
}

Member Function Documentation

void CastorClusterProducer::beginJob ( void  ) [private, virtual]

Reimplemented from edm::EDProducer.

Definition at line 241 of file CastorClusterProducer.cc.

References LogDebug.

                                     {
  LogDebug("CastorClusterProducer")
    <<"Starting CastorClusterProducer";
}
void CastorClusterProducer::endJob ( void  ) [private, virtual]

Reimplemented from edm::EDProducer.

Definition at line 247 of file CastorClusterProducer.cc.

References LogDebug.

                                   {
  LogDebug("CastorClusterProducer")
    <<"Ending CastorClusterProducer";
}
double CastorClusterProducer::phiangle ( double  testphi) [private]

Definition at line 233 of file CastorClusterProducer.cc.

References M_PI, and phi.

Referenced by produce().

                                                      {
  double phi = testphi;
  while (phi>M_PI) phi -= (2*M_PI);
  while (phi<-M_PI) phi += (2*M_PI);
  return phi;
}
void CastorClusterProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [private, virtual]

Implements edm::EDProducer.

Definition at line 108 of file CastorClusterProducer.cc.

References abs, basicjets_, reco::CastorTower::cellsBegin(), reco::CastorTower::cellsEnd(), clusteralgo_, reco::CastorTower::depth(), reco::CastorTower::emEnergy(), reco::LeafCandidate::energy(), reco::CastorTower::energy(), relval_parameters_module::energy, reco::LeafCandidate::eta(), reco::CastorTower::fhot(), edm::Event::getByLabel(), reco::Jet::getJetConstituents(), reco::CastorTower::hadEnergy(), i, input_, prof2calltree::l, LogDebug, reco::LeafCandidate::phi(), reco::CastorTower::phi(), phiangle(), position, funct::pow(), edm::RefVector< C, T, F >::push_back(), edm::Event::put(), dt_offlineAnalysis_common_cff::reco, mathSSE::sqrt(), cond::rpcobtemp::temp, and tablePrinter::width.

                                                                                 {
  
  using namespace edm;
  using namespace reco;
  using namespace TMath;
  
  LogDebug("CastorClusterProducer")
    <<"3. entering CastorClusterProducer";
  
  if ( input_ != "") {
  
  // Produce CastorClusters from CastorTowers
  
  edm::Handle<CastorTowerCollection> InputTowers;
  iEvent.getByLabel(input_, InputTowers);

  std::auto_ptr<CastorClusterCollection> OutputClustersfromClusterAlgo (new CastorClusterCollection);
  
  // get and check input size
  int nTowers = InputTowers->size();

  if (nTowers==0) LogDebug("CastorClusterProducer")<<"Warning: You are trying to run the Cluster algorithm with 0 input towers.";

  CastorTowerRefVector posInputTowers, negInputTowers;

  for (size_t i = 0; i < InputTowers->size(); ++i) {
    reco::CastorTowerRef tower_p = reco::CastorTowerRef(InputTowers, i);
    if (tower_p->eta() > 0.) posInputTowers.push_back(tower_p);
    if (tower_p->eta() < 0.) negInputTowers.push_back(tower_p);
  }
    
  // build cluster from ClusterAlgo
  if (clusteralgo_ == true) {
    // code
    iEvent.put(OutputClustersfromClusterAlgo);
  }
  
  }
  
  if ( basicjets_ != "") {
  
        Handle<BasicJetCollection> bjCollection;
        iEvent.getByLabel(basicjets_,bjCollection);
        
        Handle<CastorTowerCollection> ctCollection;
        iEvent.getByLabel("CastorTowerReco",ctCollection);
        
        std::auto_ptr<CastorClusterCollection> OutputClustersfromBasicJets (new CastorClusterCollection);
        
        if (bjCollection->size()==0) LogDebug("CastorClusterProducer")<< "Warning: You are trying to run the Cluster algorithm with 0 input basicjets.";
   
        for (unsigned i=0; i< bjCollection->size();i++) {
                const BasicJet* bj = &(*bjCollection)[i];
                
                double energy = bj->energy();
                TowerPoint temp(88.5,bj->eta(),bj->phi());
                Point position(temp);
                double emEnergy = 0.;
                double hadEnergy = 0.;
                double width = 0.;
                double depth = 0.;
                double fhot = 0.;
                double sigmaz = 0.;
                CastorTowerRefVector usedTowers;
                double zmean = 0.;
                double z2mean = 0.;
        
                std::vector<CandidatePtr> ccp = bj->getJetConstituents();
                std::vector<CandidatePtr>::const_iterator itParticle;
                for (itParticle=ccp.begin();itParticle!=ccp.end();++itParticle){            
                        const CastorTower* castorcand = dynamic_cast<const CastorTower*>(itParticle->get());
                        //cout << " castortowercandidate reference energy = " << castorcand->castorTower()->energy() << endl;
                        //cout << " castortowercandidate reference eta = " << castorcand->castorTower()->eta() << endl;
                        //cout << " castortowercandidate reference phi = " << castorcand->castorTower()->phi() << endl;
                        //cout << " castortowercandidate reference depth = " << castorcand->castorTower()->depth() << endl;
                        
                        //CastorTowerCollection *ctc = new CastorTowerCollection();
                        //ctc->push_back(*castorcand);
                        //CastorTowerRef towerref = CastorTowerRef(ctc,0);
                        
                        size_t thisone = 0;
                        for (size_t l=0;l<ctCollection->size();l++) {
                                const CastorTower ct = (*ctCollection)[l];
                                if ( std::abs(ct.phi() - castorcand->phi()) < 0.0001 ) { thisone = l;}
                        }
                        
                        CastorTowerRef towerref(ctCollection,thisone); 
                        usedTowers.push_back(towerref);
                        emEnergy += castorcand->emEnergy();
                        hadEnergy += castorcand->hadEnergy();
                        depth += castorcand->depth()*castorcand->energy();
                        width += pow(phiangle(castorcand->phi() - bj->phi()),2)*castorcand->energy();
                        fhot += castorcand->fhot()*castorcand->energy();
                        
                        // loop over cells
                        for (CastorCell_iterator it = castorcand->cellsBegin(); it != castorcand->cellsEnd(); it++) {
                                CastorCellRef cell_p = *it;
                                Point rcell = cell_p->position();
                                double Ecell = cell_p->energy();
                                zmean += Ecell*cell_p->z();
                                z2mean += Ecell*cell_p->z()*cell_p->z();
                        } // end loop over cells
                }
                //cout << "" << endl;
                
                depth = depth/energy;
                width = sqrt(width/energy);
                fhot = fhot/energy;
                
                zmean = zmean/energy;
                z2mean = z2mean/energy;
                double sigmaz2 = z2mean - zmean*zmean;
                if(sigmaz2 > 0) sigmaz = sqrt(sigmaz2);
                
                CastorCluster cc(energy,position,emEnergy,hadEnergy,emEnergy/energy,width,depth,fhot,sigmaz,usedTowers);
                OutputClustersfromBasicJets->push_back(cc);
        }
        
        iEvent.put(OutputClustersfromBasicJets);
  
  }
 
}

Member Data Documentation

std::string CastorClusterProducer::basicjets_ [private]

Definition at line 67 of file CastorClusterProducer.cc.

Referenced by produce().

Definition at line 68 of file CastorClusterProducer.cc.

Referenced by produce().

std::string CastorClusterProducer::input_ [private]

Definition at line 67 of file CastorClusterProducer.cc.

Referenced by produce().