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/CastorReco/interface/CastorCell.h"
00038 #include "DataFormats/HcalRecHit/interface/CastorRecHit.h"
00039 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00040
00041
00042
00043
00044
00045 class CastorCellProducer : public edm::EDProducer {
00046 public:
00047 explicit CastorCellProducer(const edm::ParameterSet&);
00048 ~CastorCellProducer();
00049
00050 private:
00051 virtual void beginJob() ;
00052 virtual void produce(edm::Event&, const edm::EventSetup&);
00053 virtual void endJob() ;
00054
00055
00056 typedef math::XYZPointD Point;
00057 typedef ROOT::Math::RhoZPhiPoint CellPoint;
00058 typedef std::vector<reco::CastorCell> CastorCellCollection;
00059 std::string input_;
00060 };
00061
00062
00063
00064
00065
00066 const double MYR2D = 180/M_PI;
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076 CastorCellProducer::CastorCellProducer(const edm::ParameterSet& iConfig) :
00077 input_(iConfig.getUntrackedParameter<std::string>("inputprocess","castorreco"))
00078 {
00079
00080 produces<CastorCellCollection>();
00081
00082 }
00083
00084
00085 CastorCellProducer::~CastorCellProducer()
00086 {
00087
00088
00089 }
00090
00091
00092
00093
00094
00095
00096
00097
00098 void CastorCellProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00099
00100 using namespace edm;
00101 using namespace reco;
00102 using namespace TMath;
00103
00104
00105
00106 edm::Handle<CastorRecHitCollection> InputRecHits;
00107 iEvent.getByLabel(input_, InputRecHits);
00108
00109 std::auto_ptr<CastorCellCollection> OutputCells (new CastorCellCollection);
00110
00111
00112
00113 LogDebug("CastorCellProducer")
00114 <<"1. entering CastorCellProducer ";
00115
00116 for (size_t i = 0; i < InputRecHits->size(); ++i) {
00117 const CastorRecHit & rh = (*InputRecHits)[i];
00118 int sector = rh.id().sector();
00119 int module = rh.id().module();
00120 double energy = rh.energy();
00121 int zside = rh.id().zside();
00122
00123
00124 double zCell=0.;
00125 double phiCell;
00126 double rhoCell;
00127
00128
00129 if (module < 3) {
00130
00131 if (module == 1) zCell = 14415;
00132 if (module == 2) zCell = 14464;
00133 } else {
00134
00135 zCell = 14534 + (module - 3)*92;
00136 }
00137
00138
00139 double castorphi[16];
00140 for (int j = 0; j < 16; j++) {
00141 castorphi[j] = -2.94524 + j*0.3927;
00142 }
00143 if (sector > 8) {
00144 phiCell = castorphi[sector - 9];
00145 } else {
00146 phiCell = castorphi[sector + 7];
00147 }
00148
00149
00150 if (zside <= 0) zCell = -1*zCell;
00151
00152
00153 rhoCell = 88.5;
00154
00155
00156 CellPoint tempcellposition(rhoCell,zCell,phiCell);
00157 Point cellposition(tempcellposition);
00158
00159 LogDebug("CastorCellProducer")
00160 <<"cell number: "<<i+1<<std::endl
00161 <<"rho: "<<cellposition.rho()<<" phi: "<<cellposition.phi()*MYR2D<<" eta: "<<cellposition.eta()<<std::endl
00162 <<"x: "<<cellposition.x()<<" y: "<<cellposition.y()<<" z: "<<cellposition.z();
00163
00164 if (energy > 0.) {
00165 CastorCell newCell(energy,cellposition);
00166 OutputCells->push_back(newCell);
00167 }
00168
00169 }
00170
00171 LogDebug("CastorCellProducer")
00172 <<"total number of cells in the event: "<<InputRecHits->size();
00173
00174 iEvent.put(OutputCells);
00175
00176
00177 }
00178
00179
00180 void CastorCellProducer::beginJob() {
00181 LogDebug("CastorCellProducer")
00182 <<"Starting CastorCellProducer";
00183 }
00184
00185
00186 void CastorCellProducer::endJob() {
00187 LogDebug("CastorCellProducer")
00188 <<"Ending CastorCellProducer";
00189 }
00190
00191
00192 DEFINE_FWK_MODULE(CastorCellProducer);