CMS 3D CMS Logo

CastorCellProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Castor
4 // Class: CastorCellProducer
5 //
6 
13 //
14 // Original Author: Hans Van Haevermaet, Benoit Roland
15 // Created: Wed Jul 9 14:00:40 CEST 2008
16 //
17 //
18 
19 // system include
20 #include <memory>
21 #include <vector>
22 #include <iostream>
23 #include <TMath.h>
24 #include <TRandom3.h>
25 
26 // user include
33 
34 // Castor Object include
38 
39 //
40 // class declaration
41 //
42 
44 public:
45  explicit CastorCellProducer(const edm::ParameterSet&);
46  ~CastorCellProducer() override;
47 
48 private:
49  void beginJob() override;
50  void produce(edm::Event&, const edm::EventSetup&) override;
51  void endJob() override;
52 
53  // member data
55  typedef ROOT::Math::RhoZPhiPoint CellPoint;
56  typedef std::vector<reco::CastorCell> CastorCellCollection;
58 };
59 
60 //
61 // constants, enums and typedefs
62 //
63 
64 const double MYR2D = 180 / M_PI;
65 
66 //
67 // static data member definitions
68 //
69 
70 //
71 // constructor and destructor
72 //
73 
75  tok_input_ = consumes<CastorRecHitCollection>(
76  edm::InputTag(iConfig.getUntrackedParameter<std::string>("inputprocess", "castorreco")));
77  // register your products
78  produces<CastorCellCollection>();
79  // now do what ever other initialization is needed
80 }
81 
83  // do anything here that needs to be done at desctruction time
84  // (e.g. close files, deallocate resources etc.)
85 }
86 
87 //
88 // member functions
89 //
90 
91 // ------------ method called to produce the data ------------
92 
94  using namespace edm;
95  using namespace reco;
96  using namespace TMath;
97 
98  // Produce CastorCells from CastorRecHits
99 
101  iEvent.getByToken(tok_input_, InputRecHits);
102 
103  auto OutputCells = std::make_unique<CastorCellCollection>();
104 
105  // looping over all CastorRecHits
106 
107  LogDebug("CastorCellProducer") << "1. entering CastorCellProducer ";
108 
109  for (size_t i = 0; i < InputRecHits->size(); ++i) {
110  const CastorRecHit& rh = (*InputRecHits)[i];
111  int sector = rh.id().sector();
112  int module = rh.id().module();
113  double energy = rh.energy();
114  int zside = rh.id().zside();
115 
116  // define CastorCell properties
117  double zCell = 0.;
118  double phiCell;
119  double rhoCell;
120 
121  // set z position of the cell
122  if (module < 3) {
123  // starting in EM section
124  if (module == 1)
125  zCell = 14415;
126  if (module == 2)
127  zCell = 14464;
128  } else {
129  // starting in HAD section
130  zCell = 14534 + (module - 3) * 92;
131  }
132 
133  // set phi position of the cell
134  double castorphi[16];
135  for (int j = 0; j < 16; j++) {
136  castorphi[j] = -2.94524 + j * 0.3927;
137  }
138  if (sector > 8) {
139  phiCell = castorphi[sector - 9];
140  } else {
141  phiCell = castorphi[sector + 7];
142  }
143 
144  // add condition to select in eta sides
145  if (zside <= 0)
146  zCell = -1 * zCell;
147 
148  // set rho position of the cell (inner radius 3.7cm, outer radius 14cm)
149  rhoCell = 88.5;
150 
151  // store cell position
152  CellPoint tempcellposition(rhoCell, zCell, phiCell);
153  Point cellposition(tempcellposition);
154 
155  LogDebug("CastorCellProducer") << "cell number: " << i + 1 << std::endl
156  << "rho: " << cellposition.rho() << " phi: " << cellposition.phi() * MYR2D
157  << " eta: " << cellposition.eta() << std::endl
158  << "x: " << cellposition.x() << " y: " << cellposition.y()
159  << " z: " << cellposition.z();
160 
161  if (energy > 0.) {
162  CastorCell newCell(energy, cellposition);
163  OutputCells->push_back(newCell);
164  }
165 
166  } // end loop over CastorRecHits
167 
168  LogDebug("CastorCellProducer") << "total number of cells in the event: " << InputRecHits->size();
169 
170  iEvent.put(std::move(OutputCells));
171 }
172 
173 // ------------ method called once each job just before starting event loop ------------
174 void CastorCellProducer::beginJob() { LogDebug("CastorCellProducer") << "Starting CastorCellProducer"; }
175 
176 // ------------ method called once each job just after ending the event loop ------------
177 void CastorCellProducer::endJob() { LogDebug("CastorCellProducer") << "Ending CastorCellProducer"; }
178 
179 //define this as a plug-in
mps_fire.i
i
Definition: mps_fire.py:428
CaloRecHit::energy
constexpr float energy() const
Definition: CaloRecHit.h:29
EDProducer.h
CastorCellProducer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: CastorCellProducer.cc:93
ecaldqm::zside
int zside(DetId const &)
Definition: EcalDQMCommonUtils.cc:189
edm::EDGetTokenT
Definition: EDGetToken.h:33
edm
HLT enums.
Definition: AlignableModifier.h:19
CastorCellProducer
Definition: CastorCellProducer.cc:43
CastorCell.h
CastorCellProducer::CastorCellCollection
std::vector< reco::CastorCell > CastorCellCollection
Definition: CastorCellProducer.cc:56
edm::SortedCollection::size
size_type size() const
Definition: SortedCollection.h:215
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::Handle
Definition: AssociativeIterator.h:50
CastorCellProducer::beginJob
void beginJob() override
Definition: CastorCellProducer.cc:174
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
CastorRecHit
Definition: CastorRecHit.h:7
CastorCellProducer::CastorCellProducer
CastorCellProducer(const edm::ParameterSet &)
Definition: CastorCellProducer.cc:74
CastorCellProducer::endJob
void endJob() override
Definition: CastorCellProducer.cc:177
CastorRecHit::id
HcalCastorDetId id() const
get the id
Definition: CastorRecHit.h:14
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
Point
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:57
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:223
edm::ParameterSet
Definition: ParameterSet.h:47
HcalCastorDetId::module
int module() const
get the module (1-2 for EM, 1-12 for HAD)
Definition: HcalCastorDetId.h:58
Event.h
math::XYZPointD
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double > > XYZPointD
point in space with cartesian internal representation
Definition: Point3D.h:8
HcalCastorDetId::zside
int zside() const
get the z-side of the cell (1/-1)
Definition: HcalCastorDetId.h:50
iEvent
int iEvent
Definition: GenABIO.cc:224
HcalCastorDetId::sector
int sector() const
get the sector (1-16)
Definition: HcalCastorDetId.h:62
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
edm::EventSetup
Definition: EventSetup.h:57
CastorCellProducer::~CastorCellProducer
~CastorCellProducer() override
Definition: CastorCellProducer.cc:82
eostools.move
def move(src, dest)
Definition: eostools.py:511
Frameworkfwd.h
CastorCellProducer::Point
math::XYZPointD Point
Definition: CastorCellProducer.cc:54
Point3D.h
MYR2D
const double MYR2D
Definition: CastorCellProducer.cc:64
edm::EDProducer
Definition: EDProducer.h:35
HcalRecHitCollections.h
CastorCellProducer::CellPoint
ROOT::Math::RhoZPhiPoint CellPoint
Definition: CastorCellProducer.cc:55
CastorRecHit.h
CastorCellProducer::tok_input_
edm::EDGetTokenT< CastorRecHitCollection > tok_input_
Definition: CastorCellProducer.cc:57
ParameterSet.h
reco::CastorCell
Definition: CastorCell.h:22
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15