Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #include <memory>
00023 #include <vector>
00024 #include <iostream>
00025 #include <TMath.h>
00026 #include <TRandom3.h>
00027
00028
00029 #include "FWCore/Framework/interface/Frameworkfwd.h"
00030 #include "FWCore/Framework/interface/EDProducer.h"
00031 #include "FWCore/Framework/interface/Event.h"
00032 #include "FWCore/Framework/interface/MakerMacros.h"
00033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00034 #include "DataFormats/Math/interface/Point3D.h"
00035
00036
00037 #include "DataFormats/HcalRecHit/interface/CastorRecHit.h"
00038 #include "DataFormats/HcalDetId/interface/HcalCastorDetId.h"
00039 #include "DataFormats/CastorReco/interface/CastorTower.h"
00040
00041
00042
00043
00044
00045 class CastorTowerProducer : public edm::EDProducer {
00046 public:
00047 explicit CastorTowerProducer(const edm::ParameterSet&);
00048 ~CastorTowerProducer();
00049
00050 private:
00051 virtual void beginJob() ;
00052 virtual void produce(edm::Event&, const edm::EventSetup&);
00053 virtual void endJob() ;
00054 virtual void ComputeTowerVariable(const edm::RefVector<edm::SortedCollection<CastorRecHit> >& usedRecHits, double& Ehot, double& depth);
00055
00056
00057 typedef math::XYZPointD Point;
00058 typedef ROOT::Math::RhoEtaPhiPoint TowerPoint;
00059 typedef ROOT::Math::RhoZPhiPoint CellPoint;
00060 typedef edm::SortedCollection<CastorRecHit> CastorRecHitCollection;
00061 typedef std::vector<reco::CastorTower> CastorTowerCollection;
00062 typedef edm::RefVector<CastorRecHitCollection> CastorRecHitRefVector;
00063 std::string input_;
00064 double towercut_;
00065 double mintime_;
00066 double maxtime_;
00067 };
00068
00069
00070
00071
00072
00073 const double MYR2D = 180/M_PI;
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083 CastorTowerProducer::CastorTowerProducer(const edm::ParameterSet& iConfig) :
00084 input_(iConfig.getUntrackedParameter<std::string>("inputprocess","castorreco")),
00085 towercut_(iConfig.getUntrackedParameter<double>("towercut",1.)),
00086 mintime_(iConfig.getUntrackedParameter<double>("mintime",-999)),
00087 maxtime_(iConfig.getUntrackedParameter<double>("maxtime",999))
00088 {
00089
00090 produces<CastorTowerCollection>();
00091
00092 }
00093
00094
00095 CastorTowerProducer::~CastorTowerProducer()
00096 {
00097
00098
00099 }
00100
00101
00102
00103
00104
00105
00106
00107 void CastorTowerProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00108
00109 using namespace edm;
00110 using namespace reco;
00111 using namespace TMath;
00112
00113
00114
00115 edm::Handle<CastorRecHitCollection> InputRecHits;
00116 iEvent.getByLabel(input_,InputRecHits);
00117
00118 std::auto_ptr<CastorTowerCollection> OutputTowers (new CastorTowerCollection);
00119
00120
00121 int nRecHits = InputRecHits->size();
00122
00123 LogDebug("CastorTowerProducer")
00124 <<"2. entering CastorTowerProducer"<<std::endl;
00125
00126 if (nRecHits==0)
00127 LogDebug("CastorTowerProducer") <<"Warning: You are trying to run the Tower algorithm with 0 input rechits.";
00128
00129
00130
00131
00132 double poscastortowerarray[4][16];
00133 double negcastortowerarray[4][16];
00134
00135 CastorRecHitRefVector poscastorusedrechits[16];
00136 CastorRecHitRefVector negcastorusedrechits[16];
00137
00138
00139 for (int j = 0; j < 16; j++) {
00140 poscastortowerarray[3][j] = -2.94524 + j*0.3927;
00141 poscastortowerarray[0][j] = 0.;
00142 poscastortowerarray[1][j] = 0.;
00143 poscastortowerarray[2][j] = 0.;
00144
00145 negcastortowerarray[3][j] = -2.94524 + j*0.3927;
00146 negcastortowerarray[0][j] = 0.;
00147 negcastortowerarray[1][j] = 0.;
00148 negcastortowerarray[2][j] = 0.;
00149 }
00150
00151
00152 for (unsigned int i = 0; i < InputRecHits->size(); i++) {
00153
00154 edm::Ref<CastorRecHitCollection> rechit_p = edm::Ref<CastorRecHitCollection>(InputRecHits, i);
00155
00156 double Erechit = rechit_p->energy();
00157 HcalCastorDetId id = rechit_p->id();
00158 int module = id.module();
00159 int sector = id.sector();
00160 double zrechit = 0;
00161 if (module < 3) zrechit = -14390 - 24.75 - 49.5*(module-1);
00162 if (module > 2) zrechit = -14390 - 99 - 49.5 - 99*(module-3);
00163 double phirechit = -100;
00164 if (sector < 9) phirechit = 0.19635 + (sector-1)*0.3927;
00165 if (sector > 8) phirechit = -2.94524 + (sector - 9)*0.3927;
00166
00167
00168 if (rechit_p->time() > mintime_ && rechit_p->time() < maxtime_) {
00169
00170
00171 for ( int j=0;j<16;j++) {
00172
00173
00174 if (TMath::Abs(phirechit - poscastortowerarray[3][j]) < 0.0001) {
00175
00176
00177 if (zrechit > 0.) {
00178 poscastortowerarray[0][j]+=Erechit;
00179 if (module < 3) {poscastortowerarray[1][j]+=Erechit;} else {poscastortowerarray[2][j]+=Erechit;}
00180 poscastorusedrechits[j].push_back(rechit_p);
00181 } else {
00182 negcastortowerarray[0][j]+=Erechit;
00183 if (module < 3) {negcastortowerarray[1][j]+=Erechit;} else {negcastortowerarray[2][j]+=Erechit;}
00184 negcastorusedrechits[j].push_back(rechit_p);
00185 }
00186 }
00187 }
00188 }
00189
00190 }
00191
00192
00193
00194 double fem, Ehot, depth;
00195 double rhoTower = 88.5;
00196
00197
00198 for (int k=0;k<16;k++) {
00199
00200 fem = 0;
00201 Ehot = 0;
00202 depth = 0;
00203
00204
00205 if (poscastortowerarray[0][k] > towercut_) {
00206
00207 fem = poscastortowerarray[1][k]/poscastortowerarray[0][k];
00208 CastorRecHitRefVector usedRecHits = poscastorusedrechits[k];
00209 ComputeTowerVariable(usedRecHits,Ehot,depth);
00210
00211 LogDebug("CastorTowerProducer")
00212 <<"tower "<<k+1<<": fem = "<<fem<<" ,depth = "<<depth<<" ,Ehot = "<<Ehot<<std::endl;
00213
00214 TowerPoint temptowerposition(rhoTower,5.9,poscastortowerarray[3][k]);
00215 Point towerposition(temptowerposition);
00216
00217 CastorTower newtower(poscastortowerarray[0][k],towerposition,poscastortowerarray[1][k],poscastortowerarray[2][k],fem,depth,Ehot,
00218 poscastorusedrechits[k]);
00219 OutputTowers->push_back(newtower);
00220 }
00221
00222
00223 if (negcastortowerarray[0][k] > towercut_) {
00224
00225 fem = negcastortowerarray[1][k]/negcastortowerarray[0][k];
00226 CastorRecHitRefVector usedRecHits = negcastorusedrechits[k];
00227 ComputeTowerVariable(usedRecHits,Ehot,depth);
00228
00229 LogDebug("CastorTowerProducer")
00230 <<"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;
00231
00232 TowerPoint temptowerposition(rhoTower,-5.9,negcastortowerarray[3][k]);
00233 Point towerposition(temptowerposition);
00234
00235 CastorTower newtower(negcastortowerarray[0][k],towerposition,negcastortowerarray[1][k],negcastortowerarray[2][k],fem,depth,Ehot,
00236 negcastorusedrechits[k]);
00237 OutputTowers->push_back(newtower);
00238 }
00239
00240 }
00241
00242 iEvent.put(OutputTowers);
00243 }
00244
00245
00246
00247 void CastorTowerProducer::beginJob() {
00248 LogDebug("CastorTowerProducer")
00249 <<"Starting CastorTowerProducer";
00250 }
00251
00252
00253 void CastorTowerProducer::endJob() {
00254 LogDebug("CastorTowerProducer")
00255 <<"Ending CastorTowerProducer";
00256 }
00257
00258 void CastorTowerProducer::ComputeTowerVariable(const edm::RefVector<edm::SortedCollection<CastorRecHit> >& usedRecHits, double& Ehot, double& depth) {
00259
00260 using namespace reco;
00261
00262 double Etot = 0;
00263
00264
00265 for (CastorRecHitRefVector::iterator it = usedRecHits.begin(); it != usedRecHits.end(); it++) {
00266 edm::Ref<CastorRecHitCollection> rechit_p = *it;
00267
00268 double Erechit = rechit_p->energy();
00269 HcalCastorDetId id = rechit_p->id();
00270 int module = id.module();
00271 double zrechit = 0;
00272 if (module < 3) zrechit = -14390 - 24.75 - 49.5*(module-1);
00273 if (module > 2) zrechit = -14390 - 99 - 49.5 - 99*(module-3);
00274
00275 if(Erechit > Ehot) Ehot = Erechit;
00276 depth+=Erechit*zrechit;
00277 Etot+=Erechit;
00278 }
00279
00280 depth/=Etot;
00281 Ehot/=Etot;
00282 }
00283
00284
00285 DEFINE_FWK_MODULE(CastorTowerProducer);