CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/RecoLocalCalo/Castor/src/CastorClusterProducer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    Castor
00004 // Class:      CastorClusterProducer
00005 // 
00012 //
00013 // Original Author:  Hans Van Haevermaet, Benoit Roland
00014 //         Created:  Wed Jul  9 14:00:40 CEST 2008
00015 // $Id: CastorClusterProducer.cc,v 1.9 2010/10/26 04:16:21 wmtan Exp $
00016 //
00017 //
00018 
00019 // system include 
00020 #include <memory>
00021 #include <vector>
00022 #include <iostream>
00023 #include <TMath.h>
00024 #include <TRandom3.h>
00025 
00026 // user include 
00027 #include "FWCore/Framework/interface/Frameworkfwd.h"
00028 #include "FWCore/Framework/interface/EDProducer.h"
00029 #include "FWCore/Framework/interface/Event.h"
00030 #include "FWCore/Framework/interface/MakerMacros.h"
00031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00032 #include "DataFormats/Math/interface/Point3D.h"
00033 
00034 // Castor Object include
00035 #include "DataFormats/CastorReco/interface/CastorCell.h"
00036 #include "DataFormats/CastorReco/interface/CastorTower.h"
00037 #include "DataFormats/CastorReco/interface/CastorCluster.h"
00038 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00039 #include "DataFormats/Candidate/interface/Candidate.h"
00040 #include "DataFormats/JetReco/interface/BasicJet.h"
00041 #include "DataFormats/JetReco/interface/BasicJetCollection.h"
00042 #include "DataFormats/JetReco/interface/Jet.h"
00043 
00044 
00045 
00046 //
00047 // class decleration
00048 //
00049 
00050 class CastorClusterProducer : public edm::EDProducer {
00051    public:
00052       explicit CastorClusterProducer(const edm::ParameterSet&);
00053       ~CastorClusterProducer();
00054 
00055    private:
00056       virtual void beginJob() ;
00057       virtual void produce(edm::Event&, const edm::EventSetup&);
00058       double phiangle (double testphi);
00059       virtual void endJob() ;
00060       
00061       // ----------member data ---------------------------
00062       typedef math::XYZPointD Point;
00063       typedef ROOT::Math::RhoEtaPhiPoint TowerPoint;
00064       typedef ROOT::Math::RhoZPhiPoint CellPoint;
00065       typedef std::vector<reco::CastorTower> CastorTowerCollection;
00066       typedef std::vector<reco::CastorCluster> CastorClusterCollection;
00067       std::string input_, basicjets_;
00068       bool clusteralgo_;
00069 };
00070 
00071 //
00072 // constants, enums and typedefs
00073 //
00074 
00075 const double MYR2D = 180/M_PI;
00076 
00077 //
00078 // static data member definitions
00079 //
00080 
00081 //
00082 // constructor and destructor
00083 //
00084 
00085 CastorClusterProducer::CastorClusterProducer(const edm::ParameterSet& iConfig) :
00086   input_(iConfig.getUntrackedParameter<std::string>("inputtowers","")),
00087   basicjets_(iConfig.getUntrackedParameter<std::string>("basicjets","")),
00088   clusteralgo_(iConfig.getUntrackedParameter<bool>("ClusterAlgo",false))
00089 {
00090   // register your products
00091   produces<CastorClusterCollection>();
00092   // now do what ever other initialization is needed
00093 }
00094 
00095 
00096 CastorClusterProducer::~CastorClusterProducer()
00097 {
00098   // do anything here that needs to be done at desctruction time
00099   // (e.g. close files, deallocate resources etc.)
00100 }
00101 
00102 
00103 //
00104 // member functions
00105 //
00106 
00107 // ------------ method called to produce the data  ------------
00108 void CastorClusterProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00109   
00110   using namespace edm;
00111   using namespace reco;
00112   using namespace TMath;
00113   
00114   LogDebug("CastorClusterProducer")
00115     <<"3. entering CastorClusterProducer";
00116   
00117   if ( input_ != "") {
00118   
00119   // Produce CastorClusters from CastorTowers
00120   
00121   edm::Handle<CastorTowerCollection> InputTowers;
00122   iEvent.getByLabel(input_, InputTowers);
00123 
00124   std::auto_ptr<CastorClusterCollection> OutputClustersfromClusterAlgo (new CastorClusterCollection);
00125   
00126   // get and check input size
00127   int nTowers = InputTowers->size();
00128 
00129   if (nTowers==0) LogDebug("CastorClusterProducer")<<"Warning: You are trying to run the Cluster algorithm with 0 input towers.";
00130 
00131   CastorTowerRefVector posInputTowers, negInputTowers;
00132 
00133   for (size_t i = 0; i < InputTowers->size(); ++i) {
00134     reco::CastorTowerRef tower_p = reco::CastorTowerRef(InputTowers, i);
00135     if (tower_p->eta() > 0.) posInputTowers.push_back(tower_p);
00136     if (tower_p->eta() < 0.) negInputTowers.push_back(tower_p);
00137   }
00138     
00139   // build cluster from ClusterAlgo
00140   if (clusteralgo_ == true) {
00141     // code
00142     iEvent.put(OutputClustersfromClusterAlgo);
00143   }
00144   
00145   }
00146   
00147   if ( basicjets_ != "") {
00148   
00149         Handle<BasicJetCollection> bjCollection;
00150         iEvent.getByLabel(basicjets_,bjCollection);
00151         
00152         Handle<CastorTowerCollection> ctCollection;
00153         iEvent.getByLabel("CastorTowerReco",ctCollection);
00154         
00155         std::auto_ptr<CastorClusterCollection> OutputClustersfromBasicJets (new CastorClusterCollection);
00156         
00157         if (bjCollection->size()==0) LogDebug("CastorClusterProducer")<< "Warning: You are trying to run the Cluster algorithm with 0 input basicjets.";
00158    
00159         for (unsigned i=0; i< bjCollection->size();i++) {
00160                 const BasicJet* bj = &(*bjCollection)[i];
00161                 
00162                 double energy = bj->energy();
00163                 TowerPoint temp(88.5,bj->eta(),bj->phi());
00164                 Point position(temp);
00165                 double emEnergy = 0.;
00166                 double hadEnergy = 0.;
00167                 double width = 0.;
00168                 double depth = 0.;
00169                 double fhot = 0.;
00170                 double sigmaz = 0.;
00171                 CastorTowerRefVector usedTowers;
00172                 double zmean = 0.;
00173                 double z2mean = 0.;
00174         
00175                 std::vector<CandidatePtr> ccp = bj->getJetConstituents();
00176                 std::vector<CandidatePtr>::const_iterator itParticle;
00177                 for (itParticle=ccp.begin();itParticle!=ccp.end();++itParticle){            
00178                         const CastorTower* castorcand = dynamic_cast<const CastorTower*>(itParticle->get());
00179                         //cout << " castortowercandidate reference energy = " << castorcand->castorTower()->energy() << endl;
00180                         //cout << " castortowercandidate reference eta = " << castorcand->castorTower()->eta() << endl;
00181                         //cout << " castortowercandidate reference phi = " << castorcand->castorTower()->phi() << endl;
00182                         //cout << " castortowercandidate reference depth = " << castorcand->castorTower()->depth() << endl;
00183                         
00184                         //CastorTowerCollection *ctc = new CastorTowerCollection();
00185                         //ctc->push_back(*castorcand);
00186                         //CastorTowerRef towerref = CastorTowerRef(ctc,0);
00187                         
00188                         size_t thisone = 0;
00189                         for (size_t l=0;l<ctCollection->size();l++) {
00190                                 const CastorTower ct = (*ctCollection)[l];
00191                                 if ( std::abs(ct.phi() - castorcand->phi()) < 0.0001 ) { thisone = l;}
00192                         }
00193                         
00194                         CastorTowerRef towerref(ctCollection,thisone); 
00195                         usedTowers.push_back(towerref);
00196                         emEnergy += castorcand->emEnergy();
00197                         hadEnergy += castorcand->hadEnergy();
00198                         depth += castorcand->depth()*castorcand->energy();
00199                         width += pow(phiangle(castorcand->phi() - bj->phi()),2)*castorcand->energy();
00200                         fhot += castorcand->fhot()*castorcand->energy();
00201                         
00202                         // loop over cells
00203                         for (CastorCell_iterator it = castorcand->cellsBegin(); it != castorcand->cellsEnd(); it++) {
00204                                 CastorCellRef cell_p = *it;
00205                                 Point rcell = cell_p->position();
00206                                 double Ecell = cell_p->energy();
00207                                 zmean += Ecell*cell_p->z();
00208                                 z2mean += Ecell*cell_p->z()*cell_p->z();
00209                         } // end loop over cells
00210                 }
00211                 //cout << "" << endl;
00212                 
00213                 depth = depth/energy;
00214                 width = sqrt(width/energy);
00215                 fhot = fhot/energy;
00216                 
00217                 zmean = zmean/energy;
00218                 z2mean = z2mean/energy;
00219                 double sigmaz2 = z2mean - zmean*zmean;
00220                 if(sigmaz2 > 0) sigmaz = sqrt(sigmaz2);
00221                 
00222                 CastorCluster cc(energy,position,emEnergy,hadEnergy,emEnergy/energy,width,depth,fhot,sigmaz,usedTowers);
00223                 OutputClustersfromBasicJets->push_back(cc);
00224         }
00225         
00226         iEvent.put(OutputClustersfromBasicJets);
00227   
00228   }
00229  
00230 }
00231 
00232 // help function to calculate phi within [-pi,+pi]
00233 double CastorClusterProducer::phiangle (double testphi) {
00234   double phi = testphi;
00235   while (phi>M_PI) phi -= (2*M_PI);
00236   while (phi<-M_PI) phi += (2*M_PI);
00237   return phi;
00238 }
00239 
00240 // ------------ method called once each job just before starting event loop  ------------
00241 void CastorClusterProducer::beginJob() {
00242   LogDebug("CastorClusterProducer")
00243     <<"Starting CastorClusterProducer";
00244 }
00245 
00246 // ------------ method called once each job just after ending the event loop  ------------
00247 void CastorClusterProducer::endJob() {
00248   LogDebug("CastorClusterProducer")
00249     <<"Ending CastorClusterProducer";
00250 }
00251 
00252 //define this as a plug-in
00253 DEFINE_FWK_MODULE(CastorClusterProducer);