CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
PPSPixelDigiAnalyzer Class Reference
Inheritance diagram for PPSPixelDigiAnalyzer:
edm::one::EDAnalyzer< edm::one::SharedResources > edm::one::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup) override
 
void beginJob () override
 
void endJob () override
 
 PPSPixelDigiAnalyzer (const edm::ParameterSet &pset)
 
 ~PPSPixelDigiAnalyzer () override
 
- Public Member Functions inherited from edm::one::EDAnalyzer< edm::one::SharedResources >
 EDAnalyzer ()=default
 
 EDAnalyzer (const EDAnalyzer &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
const EDAnalyzeroperator= (const EDAnalyzer &)=delete
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsInputProcessBlocks () const final
 
bool wantsProcessBlocks () const final
 
- Public Member Functions inherited from edm::one::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDAnalyzerBase () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESProxyIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex > const & esGetTokenRecordIndicesVector (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::array< std::vector< ModuleDescription const *> *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, ProductRegistry const &preg, std::map< std::string, ModuleDescription const *> const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Attributes

unsigned int cumulative_cluster_size_ [3]
 
unsigned int found_corresponding_digi_count_
 
TH2D * hAllHits
 
TH2D * hOneHitperEvent
 
TH2D * hOneHitperEvent2
 
TH2D * hOneHitperEvent2Center
 
TH2D * hOneHitperEventCenter
 
std::string label_
 
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelDigi > > pixel_token
 
edm::ESGetToken< PPSPixelTopology, PPSPixelTopologyRcdpixelTopologyToken_
 
edm::EDGetTokenT< edm::PSimHitContainerpsim_token
 
int verbosity_
 

Additional Inherited Members

- Public Types inherited from edm::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
void resetItemsToGetFrom (BranchType iType)
 

Detailed Description

Definition at line 51 of file PPSPixelDigiAnalyzer.cc.

Constructor & Destructor Documentation

◆ PPSPixelDigiAnalyzer()

PPSPixelDigiAnalyzer::PPSPixelDigiAnalyzer ( const edm::ParameterSet pset)
explicit

Definition at line 77 of file PPSPixelDigiAnalyzer.cc.

References geometryDiff::file, hAllHits, hOneHitperEvent, hOneHitperEvent2, hOneHitperEvent2Center, hOneHitperEventCenter, HLT_2022v15_cff::InputTag, label_, pixel_token, pixelTopologyToken_, muonDTDigis_cfi::pset, psim_token, and verbosity_.

78  : hAllHits(nullptr),
79  hOneHitperEvent(nullptr),
80  hOneHitperEvent2(nullptr),
81  hOneHitperEventCenter(nullptr),
82  hOneHitperEvent2Center(nullptr) {
83  label_ = pset.getUntrackedParameter<string>("label");
84  verbosity_ = pset.getParameter<int>("Verbosity");
86 #ifdef USE_MIDDLE_OF_PIXEL
87  hOneHitperEvent = file->make<TH2D>("OneHitperEvent", "One Hit per Event", 30, -8.511, -8.361, 20, 1, 1.1);
88  hOneHitperEvent2 = file->make<TH2D>("OneHitperEvent2", "One Hit per Event 2", 30, -8.511, -8.361, 20, 1, 1.1);
89 #else
90  hOneHitperEvent = file->make<TH2D>("OneHitperEvent", "One Hit per Event", 30, -8.55, -8.4, 20, 1, 1.1);
91  hOneHitperEvent2 = file->make<TH2D>("OneHitperEvent2", "One Hit per Event 2", 30, -8.55, -8.4, 20, 1, 1.1);
93  file->make<TH2D>("OneHitperEventCenter", "One Hit per Event Center", 30, -0.075, 0.075, 20, -0.05, 0.05);
95  file->make<TH2D>("OneHitperEvent2Center", "Cluster Size 2", 30, -0.075, 0.075, 20, -0.05, 0.05);
96 #endif
97  file->cd();
98  hAllHits = file->make<TH2D>("AllHits", "All Hits", 10, 1, 1.1, 10, -8.55, -8.4);
99 
100  psim_token = consumes<PSimHitContainer>(edm::InputTag("g4SimHits", "CTPPSPixelHits"));
101  pixel_token = consumes<edm::DetSetVector<CTPPSPixelDigi>>(edm::InputTag(label_, "")); //label=RPixDetDigitizer???
102  pixelTopologyToken_ = esConsumes<PPSPixelTopology, PPSPixelTopologyRcd>();
103 }
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelDigi > > pixel_token
edm::ESGetToken< PPSPixelTopology, PPSPixelTopologyRcd > pixelTopologyToken_
edm::EDGetTokenT< edm::PSimHitContainer > psim_token

◆ ~PPSPixelDigiAnalyzer()

PPSPixelDigiAnalyzer::~PPSPixelDigiAnalyzer ( )
override

Definition at line 105 of file PPSPixelDigiAnalyzer.cc.

105 {}

Member Function Documentation

◆ analyze()

void PPSPixelDigiAnalyzer::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
)
overridevirtual

Implements edm::one::EDAnalyzerBase.

Definition at line 118 of file PPSPixelDigiAnalyzer.cc.

References CENTERX, CENTERY, cumulative_cluster_size_, mps_fire::end, options_cfi::eventSetup, found_corresponding_digi_count_, hAllHits, hOneHitperEvent, hOneHitperEvent2, hOneHitperEvent2Center, hOneHitperEventCenter, pixel_token, PPSPixelTopology::pixelRange(), pixelTopologyToken_, psim_token, findQualityFiles::rr, SELECTED_PIXEL_COLUMN, SELECTED_PIXEL_ROW, SELECTED_UNITID, FastTrackerRecHitCombiner_cfi::simHits, TG184, USE_MIDDLE_OF_PIXEL_2, verbosity_, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

118  {
119  if (verbosity_ > 0)
120  edm::LogInfo("PPSPixelDigiAnalyzer") << "--- Run: " << event.id().run() << " Event: " << event.id().event();
121 
122  edm::LogInfo("PPSPixelDigiAnalyzer")
123  << " I do love Pixels ";
125  event.getByToken(psim_token, simHits);
126 
128  event.getByToken(pixel_token, CTPPSPixelDigis);
129 
130  edm::ESHandle<PPSPixelTopology> thePixelTopology = eventSetup.getHandle(pixelTopologyToken_);
131 
132  if (verbosity_ > 0)
133  edm::LogInfo("PPSPixelDigiAnalyzer") << "\n=================== RPDA Starting SimHit access"
134  << " ===================";
135 
136  if (verbosity_ > 1)
137  edm::LogInfo("PPSPixelDigiAnalyzer") << simHits->size();
138 
139  double selected_pixel_lower_x;
140  double selected_pixel_lower_y;
141  double selected_pixel_upper_x;
142  double selected_pixel_upper_y;
143  double myX = 0;
144  double myY = 0;
145 
146  thePixelTopology->pixelRange(SELECTED_PIXEL_ROW,
148  selected_pixel_lower_x,
149  selected_pixel_upper_x,
150  selected_pixel_lower_y,
151  selected_pixel_upper_y);
152 
153  double hit_inside_selected_pixel[2];
154  bool found_hit_inside_selected_pixel = false;
155 
156  for (vector<PSimHit>::const_iterator hit = simHits->begin(); hit != simHits->end(); hit++) {
157  LocalPoint entryP = hit->entryPoint();
158  LocalPoint exitP = hit->exitPoint();
159  LocalPoint midP((entryP.x() + exitP.x()) / 2., (entryP.y() + exitP.y()) / 2.);
160 
161 #ifdef USE_MIDDLE_OF_PIXEL
162  if (entryP.x() > selected_pixel_lower_x && entryP.x() < selected_pixel_upper_x &&
163  entryP.y() > (selected_pixel_lower_y + 0.115 * TG184) && entryP.y() < (selected_pixel_upper_y + 0.115 * TG184)
164 #else
166  if (midP.x() > selected_pixel_lower_x && midP.x() < selected_pixel_upper_x && midP.y() > selected_pixel_lower_y &&
167  midP.y() < selected_pixel_upper_y
168 #else
169  if (entryP.x() > selected_pixel_lower_x && entryP.x() < selected_pixel_upper_x &&
170  entryP.y() > selected_pixel_lower_y && entryP.y() < selected_pixel_upper_y
171 #endif
172 #endif
173  && hit->detUnitId() == SELECTED_UNITID) {
174  hit_inside_selected_pixel[0] = entryP.x();
175  hit_inside_selected_pixel[1] = entryP.y();
176  found_hit_inside_selected_pixel = true;
177 #ifdef USE_MIDDLE_OF_PIXEL_2
178  hAllHits->Fill(midP.x(), midP.y());
179  myX = midP.x();
180  myY = midP.y();
181 #else
182  hAllHits->Fill(entryP.x(), entryP.y());
183  myX = entryP.x();
184  myY = entryP.y();
185 #endif
186  if (verbosity_ > 2)
187  edm::LogInfo("PPSPixelDigiAnalyzer") << hit_inside_selected_pixel[0] << " " << hit_inside_selected_pixel[1];
188  }
189 
190  //--------------
191 
192  if (verbosity_ > 1)
193  if (hit->timeOfFlight() > 0) {
194  edm::LogInfo("PPSPixelDigiAnalyzer")
195  << "DetId: " << hit->detUnitId() << "PID: " << hit->particleType() << " TOF: " << hit->timeOfFlight()
196  << " Proc Type: " << hit->processType() << " p: " << hit->pabs() << " x = " << entryP.x()
197  << " y = " << entryP.y() << " z = " << entryP.z();
198  }
199  }
200 
201  if (verbosity_ > 0)
202  edm::LogInfo("PPSPixelDigiAnalyzer") << "\n=================== RPDA Starting Digi access"
203  << " ===================";
204  int numberOfDetUnits = 0;
205 
206  // Iterate on detector units
207  edm::DetSetVector<CTPPSPixelDigi>::const_iterator DSViter = CTPPSPixelDigis->begin();
208 
209  for (; DSViter != CTPPSPixelDigis->end(); DSViter++) {
210  ++numberOfDetUnits;
211 
212  DetId detIdObject(DSViter->detId());
213  if (verbosity_ > 1)
214  edm::LogInfo("PPSPixelDigiAnalyzer") << "DetId: " << DSViter->detId();
215 
216  bool found_corresponding_digi = false;
217  unsigned int corresponding_digi_cluster_size = 0;
218 
219  // looping over digis in a unit id
222 
223  if (verbosity_ > 2) {
224  edm::LogInfo("PPSPixelDigiAnalyzer") << "FF " << DSViter->detId();
225  for (edm::DetSet<CTPPSPixelDigi>::const_iterator di = begin; di != end; di++) {
226  edm::LogInfo("PPSPixelDigiAnalyzer") << " Digi row " << di->row() << ", col " << di->column();
227 
228  // reconvert the digi to local coordinates
229  double lx;
230  double ly;
231  double ux;
232  double uy;
233  unsigned int rr = di->row();
234  unsigned int cc = di->column();
235  thePixelTopology->pixelRange(rr, cc, lx, ux, ly, uy);
236 
237  edm::LogInfo("PPSPixelDigiAnalyzer")
238  << " pixel boundaries x low up, y low up " << lx << " " << ux << " " << ly << " " << uy;
239  }
240  }
241  if (DSViter->detId() == SELECTED_UNITID && found_hit_inside_selected_pixel) {
242  for (edm::DetSet<CTPPSPixelDigi>::const_iterator di = begin; di != end; di++) {
243  if (verbosity_ > 1)
244  edm::LogInfo("PPSPixelDigiAnalyzer") << " Digi row " << di->row() << ", col " << di->column();
245 
246  if (di->row() == SELECTED_PIXEL_ROW && di->column() == SELECTED_PIXEL_COLUMN) {
248  found_corresponding_digi = true;
249  corresponding_digi_cluster_size = 1;
250  }
251  }
252  //if coresponding digi found, re-loop to look for adjacent pixels
253  if (found_corresponding_digi) {
254  for (edm::DetSet<CTPPSPixelDigi>::const_iterator di = begin; di != end; di++) {
255  if (verbosity_ > 1)
256  edm::LogInfo("PPSPixelDigiAnalyzer") << " Digi row " << di->row() << ", col " << di->column();
257 
258  if ((di->row() == SELECTED_PIXEL_ROW + 1 && di->column() == SELECTED_PIXEL_COLUMN) ||
259  (di->row() == SELECTED_PIXEL_ROW - 1 && di->column() == SELECTED_PIXEL_COLUMN) ||
260  (di->row() == SELECTED_PIXEL_ROW && di->column() == SELECTED_PIXEL_COLUMN + 1) ||
261  (di->row() == SELECTED_PIXEL_ROW && di->column() == SELECTED_PIXEL_COLUMN - 1)) {
262  corresponding_digi_cluster_size++;
263  edm::LogInfo("PPSPixelDigiAnalyzer") << " Digi row " << di->row() << ", col " << di->column();
264  }
265  }
266  }
267  }
268  if (corresponding_digi_cluster_size > 0) {
269  edm::LogInfo("PPSPixelDigiAnalyzer")
270  << "corresponding_digi_cluster_size in the event: " << corresponding_digi_cluster_size;
271  hOneHitperEvent->Fill(myY, myX);
272  hOneHitperEventCenter->Fill(myY - CENTERY, myX - CENTERX);
273  if (corresponding_digi_cluster_size < 3) {
274  cumulative_cluster_size_[corresponding_digi_cluster_size - 1]++;
275  if (corresponding_digi_cluster_size > 1) {
276  hOneHitperEvent2->Fill(myY, myX);
277  hOneHitperEvent2Center->Fill(myY - CENTERY, myX - CENTERX);
278  }
279  } else {
281  hOneHitperEvent2->Fill(myY, myX);
282  hOneHitperEvent2Center->Fill(myY - CENTERY, myX - CENTERX);
283  }
284  }
285  }
286 
287  if (verbosity_ > 1)
288  edm::LogInfo("PPSPixelDigiAnalyzer") << "numberOfDetUnits in the event: " << numberOfDetUnits;
289 }
unsigned int cumulative_cluster_size_[3]
T z() const
Definition: PV3DBase.h:61
unsigned int found_corresponding_digi_count_
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
#define SELECTED_PIXEL_ROW
#define TG184
#define CENTERX
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelDigi > > pixel_token
#define CENTERY
edm::ESGetToken< PPSPixelTopology, PPSPixelTopologyRcd > pixelTopologyToken_
void pixelRange(unsigned int arow, unsigned int acol, double &lower_x, double &higher_x, double &lower_y, double &higher_y) const
Log< level::Info, false > LogInfo
Definition: DetId.h:17
#define SELECTED_PIXEL_COLUMN
edm::EDGetTokenT< edm::PSimHitContainer > psim_token
#define SELECTED_UNITID
#define USE_MIDDLE_OF_PIXEL_2
collection_type::const_iterator const_iterator
Definition: DetSet.h:31
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:102

◆ beginJob()

void PPSPixelDigiAnalyzer::beginJob ( void  )
overridevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 107 of file PPSPixelDigiAnalyzer.cc.

References a, cumulative_cluster_size_, and found_corresponding_digi_count_.

107  {
109  for (int a = 0; a < 3; a++)
111 }
unsigned int cumulative_cluster_size_[3]
unsigned int found_corresponding_digi_count_
double a
Definition: hdecay.h:119

◆ endJob()

void PPSPixelDigiAnalyzer::endJob ( void  )
overridevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 112 of file PPSPixelDigiAnalyzer.cc.

References cumulative_cluster_size_, and found_corresponding_digi_count_.

112  {
113  edm::LogInfo("PPSPixelDigiAnalyzer") << "found_corresponding_digi_count_: " << found_corresponding_digi_count_;
114  edm::LogInfo("PPSPixelDigiAnalyzer") << "Cumulative cluster size (1,2,>2) = " << cumulative_cluster_size_[0] << ", "
116 }
unsigned int cumulative_cluster_size_[3]
unsigned int found_corresponding_digi_count_
Log< level::Info, false > LogInfo

Member Data Documentation

◆ cumulative_cluster_size_

unsigned int PPSPixelDigiAnalyzer::cumulative_cluster_size_[3]
private

Definition at line 74 of file PPSPixelDigiAnalyzer.cc.

Referenced by analyze(), beginJob(), and endJob().

◆ found_corresponding_digi_count_

unsigned int PPSPixelDigiAnalyzer::found_corresponding_digi_count_
private

Definition at line 73 of file PPSPixelDigiAnalyzer.cc.

Referenced by analyze(), beginJob(), and endJob().

◆ hAllHits

TH2D* PPSPixelDigiAnalyzer::hAllHits
private

Definition at line 60 of file PPSPixelDigiAnalyzer.cc.

Referenced by analyze(), and PPSPixelDigiAnalyzer().

◆ hOneHitperEvent

TH2D* PPSPixelDigiAnalyzer::hOneHitperEvent
private

Definition at line 61 of file PPSPixelDigiAnalyzer.cc.

Referenced by analyze(), and PPSPixelDigiAnalyzer().

◆ hOneHitperEvent2

TH2D* PPSPixelDigiAnalyzer::hOneHitperEvent2
private

Definition at line 62 of file PPSPixelDigiAnalyzer.cc.

Referenced by analyze(), and PPSPixelDigiAnalyzer().

◆ hOneHitperEvent2Center

TH2D* PPSPixelDigiAnalyzer::hOneHitperEvent2Center
private

Definition at line 64 of file PPSPixelDigiAnalyzer.cc.

Referenced by analyze(), and PPSPixelDigiAnalyzer().

◆ hOneHitperEventCenter

TH2D* PPSPixelDigiAnalyzer::hOneHitperEventCenter
private

Definition at line 63 of file PPSPixelDigiAnalyzer.cc.

Referenced by analyze(), and PPSPixelDigiAnalyzer().

◆ label_

std::string PPSPixelDigiAnalyzer::label_
private

◆ pixel_token

edm::EDGetTokenT<edm::DetSetVector<CTPPSPixelDigi> > PPSPixelDigiAnalyzer::pixel_token
private

Definition at line 70 of file PPSPixelDigiAnalyzer.cc.

Referenced by analyze(), and PPSPixelDigiAnalyzer().

◆ pixelTopologyToken_

edm::ESGetToken<PPSPixelTopology, PPSPixelTopologyRcd> PPSPixelDigiAnalyzer::pixelTopologyToken_
private

Definition at line 71 of file PPSPixelDigiAnalyzer.cc.

Referenced by analyze(), and PPSPixelDigiAnalyzer().

◆ psim_token

edm::EDGetTokenT<edm::PSimHitContainer> PPSPixelDigiAnalyzer::psim_token
private

Definition at line 69 of file PPSPixelDigiAnalyzer.cc.

Referenced by analyze(), and PPSPixelDigiAnalyzer().

◆ verbosity_

int PPSPixelDigiAnalyzer::verbosity_
private

Definition at line 68 of file PPSPixelDigiAnalyzer.cc.

Referenced by analyze(), and PPSPixelDigiAnalyzer().