CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/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.10 2011/02/24 09:36:46 hvanhaev 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/HcalRecHit/interface/CastorRecHit.h"
00036 #include "DataFormats/HcalDetId/interface/HcalCastorDetId.h"
00037 #include "DataFormats/CastorReco/interface/CastorTower.h"
00038 #include "DataFormats/CastorReco/interface/CastorCluster.h"
00039 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00040 #include "DataFormats/Candidate/interface/Candidate.h"
00041 #include "DataFormats/JetReco/interface/BasicJet.h"
00042 #include "DataFormats/JetReco/interface/BasicJetCollection.h"
00043 #include "DataFormats/JetReco/interface/Jet.h"
00044 
00045 
00046 
00047 //
00048 // class decleration
00049 //
00050 
00051 class CastorClusterProducer : public edm::EDProducer {
00052    public:
00053       explicit CastorClusterProducer(const edm::ParameterSet&);
00054       ~CastorClusterProducer();
00055 
00056    private:
00057       virtual void beginJob() ;
00058       virtual void produce(edm::Event&, const edm::EventSetup&);
00059       double phiangle (double testphi);
00060       virtual void endJob() ;
00061       
00062       // ----------member data ---------------------------
00063       typedef math::XYZPointD Point;
00064       typedef ROOT::Math::RhoEtaPhiPoint TowerPoint;
00065       typedef ROOT::Math::RhoZPhiPoint CellPoint;
00066       typedef std::vector<reco::CastorTower> CastorTowerCollection;
00067       typedef std::vector<reco::CastorCluster> CastorClusterCollection;
00068       std::string input_, basicjets_;
00069       bool clusteralgo_;
00070 };
00071 
00072 //
00073 // constants, enums and typedefs
00074 //
00075 
00076 const double MYR2D = 180/M_PI;
00077 
00078 //
00079 // static data member definitions
00080 //
00081 
00082 //
00083 // constructor and destructor
00084 //
00085 
00086 CastorClusterProducer::CastorClusterProducer(const edm::ParameterSet& iConfig) :
00087   input_(iConfig.getUntrackedParameter<std::string>("inputtowers","")),
00088   basicjets_(iConfig.getUntrackedParameter<std::string>("basicjets","")),
00089   clusteralgo_(iConfig.getUntrackedParameter<bool>("ClusterAlgo",false))
00090 {
00091   // register your products
00092   produces<CastorClusterCollection>();
00093   // now do what ever other initialization is needed
00094 }
00095 
00096 
00097 CastorClusterProducer::~CastorClusterProducer()
00098 {
00099   // do anything here that needs to be done at desctruction time
00100   // (e.g. close files, deallocate resources etc.)
00101 }
00102 
00103 
00104 //
00105 // member functions
00106 //
00107 
00108 // ------------ method called to produce the data  ------------
00109 void CastorClusterProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00110   
00111   using namespace edm;
00112   using namespace reco;
00113   using namespace TMath;
00114   
00115   LogDebug("CastorClusterProducer")
00116     <<"3. entering CastorClusterProducer";
00117   
00118   if ( input_ != "") {
00119   
00120   // Produce CastorClusters from CastorTowers
00121   
00122   edm::Handle<CastorTowerCollection> InputTowers;
00123   iEvent.getByLabel(input_, InputTowers);
00124 
00125   std::auto_ptr<CastorClusterCollection> OutputClustersfromClusterAlgo (new CastorClusterCollection);
00126   
00127   // get and check input size
00128   int nTowers = InputTowers->size();
00129 
00130   if (nTowers==0) LogDebug("CastorClusterProducer")<<"Warning: You are trying to run the Cluster algorithm with 0 input towers.";
00131 
00132   CastorTowerRefVector posInputTowers, negInputTowers;
00133 
00134   for (size_t i = 0; i < InputTowers->size(); ++i) {
00135     reco::CastorTowerRef tower_p = reco::CastorTowerRef(InputTowers, i);
00136     if (tower_p->eta() > 0.) posInputTowers.push_back(tower_p);
00137     if (tower_p->eta() < 0.) negInputTowers.push_back(tower_p);
00138   }
00139     
00140   // build cluster from ClusterAlgo
00141   if (clusteralgo_ == true) {
00142     // code
00143     iEvent.put(OutputClustersfromClusterAlgo);
00144   }
00145   
00146   }
00147   
00148   if ( basicjets_ != "") {
00149   
00150         Handle<BasicJetCollection> bjCollection;
00151         iEvent.getByLabel(basicjets_,bjCollection);
00152         
00153         Handle<CastorTowerCollection> ctCollection;
00154         iEvent.getByLabel("CastorTowerReco",ctCollection);
00155         
00156         std::auto_ptr<CastorClusterCollection> OutputClustersfromBasicJets (new CastorClusterCollection);
00157         
00158         if (bjCollection->size()==0) LogDebug("CastorClusterProducer")<< "Warning: You are trying to run the Cluster algorithm with 0 input basicjets.";
00159    
00160         for (unsigned i=0; i< bjCollection->size();i++) {
00161                 const BasicJet* bj = &(*bjCollection)[i];
00162                 
00163                 double energy = bj->energy();
00164                 TowerPoint temp(88.5,bj->eta(),bj->phi());
00165                 Point position(temp);
00166                 double emEnergy = 0.;
00167                 double hadEnergy = 0.;
00168                 double width = 0.;
00169                 double depth = 0.;
00170                 double fhot = 0.;
00171                 double sigmaz = 0.;
00172                 CastorTowerRefVector usedTowers;
00173                 double zmean = 0.;
00174                 double z2mean = 0.;
00175         
00176                 std::vector<CandidatePtr> ccp = bj->getJetConstituents();
00177                 std::vector<CandidatePtr>::const_iterator itParticle;
00178                 for (itParticle=ccp.begin();itParticle!=ccp.end();++itParticle){            
00179                         const CastorTower* castorcand = dynamic_cast<const CastorTower*>(itParticle->get());
00180                         //cout << " castortowercandidate reference energy = " << castorcand->castorTower()->energy() << endl;
00181                         //cout << " castortowercandidate reference eta = " << castorcand->castorTower()->eta() << endl;
00182                         //cout << " castortowercandidate reference phi = " << castorcand->castorTower()->phi() << endl;
00183                         //cout << " castortowercandidate reference depth = " << castorcand->castorTower()->depth() << endl;
00184                         
00185                         //CastorTowerCollection *ctc = new CastorTowerCollection();
00186                         //ctc->push_back(*castorcand);
00187                         //CastorTowerRef towerref = CastorTowerRef(ctc,0);
00188                         
00189                         size_t thisone = 0;
00190                         for (size_t l=0;l<ctCollection->size();l++) {
00191                                 const CastorTower ct = (*ctCollection)[l];
00192                                 if ( std::abs(ct.phi() - castorcand->phi()) < 0.0001 ) { thisone = l;}
00193                         }
00194                         
00195                         CastorTowerRef towerref(ctCollection,thisone); 
00196                         usedTowers.push_back(towerref);
00197                         emEnergy += castorcand->emEnergy();
00198                         hadEnergy += castorcand->hadEnergy();
00199                         depth += castorcand->depth()*castorcand->energy();
00200                         width += pow(phiangle(castorcand->phi() - bj->phi()),2)*castorcand->energy();
00201                         fhot += castorcand->fhot()*castorcand->energy();
00202                         
00203                         // loop over rechits
00204                         for (edm::RefVector<edm::SortedCollection<CastorRecHit> >::iterator it = castorcand->rechitsBegin(); it != castorcand->rechitsEnd(); it++) {
00205                                  edm::Ref<edm::SortedCollection<CastorRecHit> > rechit_p = *it;                         
00206                                  double Erechit = rechit_p->energy();
00207                                  HcalCastorDetId id = rechit_p->id();
00208                                  int module = id.module();                                      
00209                                  double zrechit = 0;     
00210                                  if (module < 3) zrechit = -14390 - 24.75 - 49.5*(module-1);     
00211                                  if (module > 2) zrechit = -14390 - 99 - 49.5 - 99*(module-3);   
00212                                  zmean += Erechit*zrechit;       
00213                                  z2mean += Erechit*zrechit*zrechit;
00214                         } // end loop over rechits
00215                 }
00216                 //cout << "" << endl;
00217                 
00218                 depth = depth/energy;
00219                 width = sqrt(width/energy);
00220                 fhot = fhot/energy;
00221                 
00222                 zmean = zmean/energy;
00223                 z2mean = z2mean/energy;
00224                 double sigmaz2 = z2mean - zmean*zmean;
00225                 if(sigmaz2 > 0) sigmaz = sqrt(sigmaz2);
00226                 
00227                 CastorCluster cc(energy,position,emEnergy,hadEnergy,emEnergy/energy,width,depth,fhot,sigmaz,usedTowers);
00228                 OutputClustersfromBasicJets->push_back(cc);
00229         }
00230         
00231         iEvent.put(OutputClustersfromBasicJets);
00232   
00233   }
00234  
00235 }
00236 
00237 // help function to calculate phi within [-pi,+pi]
00238 double CastorClusterProducer::phiangle (double testphi) {
00239   double phi = testphi;
00240   while (phi>M_PI) phi -= (2*M_PI);
00241   while (phi<-M_PI) phi += (2*M_PI);
00242   return phi;
00243 }
00244 
00245 // ------------ method called once each job just before starting event loop  ------------
00246 void CastorClusterProducer::beginJob() {
00247   LogDebug("CastorClusterProducer")
00248     <<"Starting CastorClusterProducer";
00249 }
00250 
00251 // ------------ method called once each job just after ending the event loop  ------------
00252 void CastorClusterProducer::endJob() {
00253   LogDebug("CastorClusterProducer")
00254     <<"Ending CastorClusterProducer";
00255 }
00256 
00257 //define this as a plug-in
00258 DEFINE_FWK_MODULE(CastorClusterProducer);