CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CastorClusterProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Castor
4 // Class: CastorClusterProducer
5 //
12 //
13 // Original Author: Hans Van Haevermaet, Benoit Roland
14 // Created: Wed Jul 9 14:00:40 CEST 2008
15 //
16 //
17 
18 // system include
19 #include <memory>
20 #include <vector>
21 #include <iostream>
22 #include <TMath.h>
23 #include <TRandom3.h>
24 
25 // user include
32 
33 // Castor Object include
43 
44 
45 
46 //
47 // class decleration
48 //
49 
51  public:
54 
55  private:
56  virtual void beginJob() override ;
57  virtual void produce(edm::Event&, const edm::EventSetup&) override;
58  double phiangle (double testphi);
59  virtual void endJob() override ;
60 
61  // ----------member data ---------------------------
63  typedef ROOT::Math::RhoEtaPhiPoint TowerPoint;
64  typedef ROOT::Math::RhoZPhiPoint CellPoint;
65  typedef std::vector<reco::CastorTower> CastorTowerCollection;
66  typedef std::vector<reco::CastorCluster> CastorClusterCollection;
72 };
73 
74 //
75 // constants, enums and typedefs
76 //
77 
78 const double MYR2D = 180/M_PI;
79 
80 //
81 // static data member definitions
82 //
83 
84 //
85 // constructor and destructor
86 //
87 
89  input_(iConfig.getUntrackedParameter<std::string>("inputtowers","")),
90  basicjets_(iConfig.getUntrackedParameter<std::string>("basicjets","")),
91  clusteralgo_(iConfig.getUntrackedParameter<bool>("ClusterAlgo",false))
92 {
93  // register for data access
94  tok_input_ = consumes<CastorTowerCollection>(edm::InputTag(input_));
95  tok_jets_ = consumes<reco::BasicJetCollection>(edm::InputTag(basicjets_));
96  tok_tower_ = consumes<CastorTowerCollection>(edm::InputTag("CastorTowerReco"));
97  // register your products
98  produces<CastorClusterCollection>();
99  // now do what ever other initialization is needed
100 }
101 
102 
104 {
105  // do anything here that needs to be done at desctruction time
106  // (e.g. close files, deallocate resources etc.)
107 }
108 
109 
110 //
111 // member functions
112 //
113 
114 // ------------ method called to produce the data ------------
116 
117  using namespace edm;
118  using namespace reco;
119  using namespace TMath;
120 
121  LogDebug("CastorClusterProducer")
122  <<"3. entering CastorClusterProducer";
123 
124  if ( input_ != "") {
125 
126  // Produce CastorClusters from CastorTowers
127 
129  iEvent.getByToken(tok_input_, InputTowers);
130 
131  std::auto_ptr<CastorClusterCollection> OutputClustersfromClusterAlgo (new CastorClusterCollection);
132 
133  // get and check input size
134  int nTowers = InputTowers->size();
135 
136  if (nTowers==0) LogDebug("CastorClusterProducer")<<"Warning: You are trying to run the Cluster algorithm with 0 input towers.";
137 
138  CastorTowerRefVector posInputTowers, negInputTowers;
139 
140  for (size_t i = 0; i < InputTowers->size(); ++i) {
141  reco::CastorTowerRef tower_p = reco::CastorTowerRef(InputTowers, i);
142  if (tower_p->eta() > 0.) posInputTowers.push_back(tower_p);
143  if (tower_p->eta() < 0.) negInputTowers.push_back(tower_p);
144  }
145 
146  // build cluster from ClusterAlgo
147  if (clusteralgo_ == true) {
148  // code
149  iEvent.put(OutputClustersfromClusterAlgo);
150  }
151 
152  }
153 
154  if ( basicjets_ != "") {
155 
156  Handle<BasicJetCollection> bjCollection;
157  iEvent.getByToken(tok_jets_,bjCollection);
158 
159  Handle<CastorTowerCollection> ctCollection;
160  iEvent.getByToken(tok_tower_,ctCollection);
161 
162  std::auto_ptr<CastorClusterCollection> OutputClustersfromBasicJets (new CastorClusterCollection);
163 
164  if (bjCollection->size()==0) LogDebug("CastorClusterProducer")<< "Warning: You are trying to run the Cluster algorithm with 0 input basicjets.";
165 
166  for (unsigned i=0; i< bjCollection->size();i++) {
167  const BasicJet* bj = &(*bjCollection)[i];
168 
169  double energy = bj->energy();
170  TowerPoint temp(88.5,bj->eta(),bj->phi());
172  double emEnergy = 0.;
173  double hadEnergy = 0.;
174  double width = 0.;
175  double depth = 0.;
176  double fhot = 0.;
177  double sigmaz = 0.;
178  CastorTowerRefVector usedTowers;
179  double zmean = 0.;
180  double z2mean = 0.;
181 
182  std::vector<CandidatePtr> ccp = bj->getJetConstituents();
183  std::vector<CandidatePtr>::const_iterator itParticle;
184  for (itParticle=ccp.begin();itParticle!=ccp.end();++itParticle){
185  const CastorTower* castorcand = dynamic_cast<const CastorTower*>(itParticle->get());
186  //cout << " castortowercandidate reference energy = " << castorcand->castorTower()->energy() << endl;
187  //cout << " castortowercandidate reference eta = " << castorcand->castorTower()->eta() << endl;
188  //cout << " castortowercandidate reference phi = " << castorcand->castorTower()->phi() << endl;
189  //cout << " castortowercandidate reference depth = " << castorcand->castorTower()->depth() << endl;
190 
191  //CastorTowerCollection *ctc = new CastorTowerCollection();
192  //ctc->push_back(*castorcand);
193  //CastorTowerRef towerref = CastorTowerRef(ctc,0);
194 
195  size_t thisone = 0;
196  for (size_t l=0;l<ctCollection->size();l++) {
197  const CastorTower ct = (*ctCollection)[l];
198  if ( std::abs(ct.phi() - castorcand->phi()) < 0.0001 ) { thisone = l;}
199  }
200 
201  CastorTowerRef towerref(ctCollection,thisone);
202  usedTowers.push_back(towerref);
203  emEnergy += castorcand->emEnergy();
204  hadEnergy += castorcand->hadEnergy();
205  depth += castorcand->depth()*castorcand->energy();
206  width += pow(phiangle(castorcand->phi() - bj->phi()),2)*castorcand->energy();
207  fhot += castorcand->fhot()*castorcand->energy();
208 
209  // loop over rechits
210  for (edm::RefVector<edm::SortedCollection<CastorRecHit> >::iterator it = castorcand->rechitsBegin(); it != castorcand->rechitsEnd(); it++) {
212  double Erechit = rechit_p->energy();
213  HcalCastorDetId id = rechit_p->id();
214  int module = id.module();
215  double zrechit = 0;
216  if (module < 3) zrechit = -14390 - 24.75 - 49.5*(module-1);
217  if (module > 2) zrechit = -14390 - 99 - 49.5 - 99*(module-3);
218  zmean += Erechit*zrechit;
219  z2mean += Erechit*zrechit*zrechit;
220  } // end loop over rechits
221  }
222  //cout << "" << endl;
223 
224  depth = depth/energy;
225  width = sqrt(width/energy);
226  fhot = fhot/energy;
227 
228  zmean = zmean/energy;
229  z2mean = z2mean/energy;
230  double sigmaz2 = z2mean - zmean*zmean;
231  if(sigmaz2 > 0) sigmaz = sqrt(sigmaz2);
232 
233  CastorCluster cc(energy,position,emEnergy,hadEnergy,emEnergy/energy,width,depth,fhot,sigmaz,usedTowers);
234  OutputClustersfromBasicJets->push_back(cc);
235  }
236 
237  iEvent.put(OutputClustersfromBasicJets);
238 
239  }
240 
241 }
242 
243 // help function to calculate phi within [-pi,+pi]
244 double CastorClusterProducer::phiangle (double testphi) {
245  double phi = testphi;
246  while (phi>M_PI) phi -= (2*M_PI);
247  while (phi<-M_PI) phi += (2*M_PI);
248  return phi;
249 }
250 
251 // ------------ method called once each job just before starting event loop ------------
253  LogDebug("CastorClusterProducer")
254  <<"Starting CastorClusterProducer";
255 }
256 
257 // ------------ method called once each job just after ending the event loop ------------
259  LogDebug("CastorClusterProducer")
260  <<"Ending CastorClusterProducer";
261 }
262 
263 //define this as a plug-in
#define LogDebug(id)
edm::EDGetTokenT< CastorTowerCollection > tok_input_
double emEnergy() const
tower em energy
Definition: CastorTower.h:48
virtual double energy() const GCC11_FINAL
energy
int i
Definition: DBlmapReader.cc:9
virtual void beginJob() override
module()
Definition: vlib.cc:994
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< CastorTowerCollection > tok_tower_
virtual void endJob() override
double phiangle(double testphi)
edm::Ref< CastorTowerCollection > CastorTowerRef
Definition: CastorTower.h:129
CastorRecHitRefs::iterator rechitsBegin() const
fist iterator over CastorRecHit constituents
Definition: CastorTower.h:66
Jets made from CaloTowers.
Definition: BasicJet.h:20
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
CastorClusterProducer(const edm::ParameterSet &)
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
std::vector< reco::CastorCluster > CastorClusterCollection
T sqrt(T t)
Definition: SSEVec.h:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< reco::CastorTower > CastorTowerCollection
edm::EDGetTokenT< reco::BasicJetCollection > tok_jets_
ROOT::Math::RhoEtaPhiPoint TowerPoint
virtual void produce(edm::Event &, const edm::EventSetup &) override
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
#define M_PI
Definition: BFit3D.cc:3
ROOT::Math::RhoZPhiPoint CellPoint
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double > > XYZPointD
point in space with cartesian internal representation
Definition: Point3D.h:8
const double MYR2D
double fhot() const
tower hotcell/tot ratio
Definition: CastorTower.h:60
double depth() const
tower depth in z
Definition: CastorTower.h:57
double hadEnergy() const
tower had energy
Definition: CastorTower.h:51
CastorRecHitRefs::iterator rechitsEnd() const
last iterator over CastorRecHit constituents
Definition: CastorTower.h:69
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:64
ProductID id() const
Accessor for product ID.
Definition: Ref.h:256
volatile std::atomic< bool > shutdown_flag false
Definition: vlib.h:208
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
virtual Constituents getJetConstituents() const
list of constituents
Definition: Jet.cc:350
Definition: DDAxes.h:10