CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes
TTStubBuilder< T > Class Template Reference

Plugin to load the Stub finding algorithm and produce the collection of Stubs that goes in the event content. More...

#include <TTStubBuilder.h>

Inheritance diagram for TTStubBuilder< T >:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 TTStubBuilder (const edm::ParameterSet &iConfig)
 Constructor. More...
 
 ~TTStubBuilder () override
 Destructor;. More...
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDProducer () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () noexcept(false) 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
 
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::vector< ModuleDescription const * > &modules, 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
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void beginRun (const edm::Run &run, const edm::EventSetup &iSetup) override
 Mandatory methods. More...
 
void endRun (const edm::Run &run, const edm::EventSetup &iSetup) override
 End run. More...
 
template<>
void produce (edm::Event &iEvent, const edm::EventSetup &iSetup)
 Implementation of methods of TTClusterBuilder.h. More...
 
void produce (edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
template<>
void produce (edm::Event &iEvent, const edm::EventSetup &iSetup)
 Implement the producer. More...
 

Static Private Member Functions

static bool SortStubBendPairs (const std::pair< unsigned int, double > &left, const std::pair< unsigned int, double > &right)
 Sort routine for stub ordering. More...
 
static bool SortStubsBend (const TTStub< T > &left, const TTStub< T > &right)
 Analogous sorting routine directly from stubs. More...
 

Private Attributes

bool applyFE
 
edm::EDGetTokenT< edmNew::DetSetVector< TTCluster< T > > > clustersToken
 
bool ForbidMultipleStubs
 
unsigned int high_rate_max_ring [5]
 
int ievt
 
unsigned int maxStubs_2S
 
unsigned int maxStubs_2S_CIC_5
 
unsigned int maxStubs_PS
 
unsigned int maxStubs_PS_CIC_10
 
unsigned int maxStubs_PS_CIC_5
 
std::unordered_map< int, int > moduleStubs_CBC
 
std::unordered_map< int, std::vector< TTStub< Ref_Phase2TrackerDigi_ > > > moduleStubs_CIC
 Temporary storage for stubs before max check. More...
 
std::unordered_map< int, int > moduleStubs_MPA
 
unsigned int tedd1_maxring
 
unsigned int tedd2_maxring
 
edm::ESHandle< TTStubAlgorithm< T > > theStubFindingAlgoHandle
 Data members. More...
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::ProducerBase
ProducesCollector producesCollector ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
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<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)
 

Detailed Description

template<typename T>
class TTStubBuilder< T >

Plugin to load the Stub finding algorithm and produce the collection of Stubs that goes in the event content.

After moving from SimDataFormats to DataFormats, the template structure of the class was maintained in order to accomodate any types other than PixelDigis in case there is such a need in the future.

Author
Andrew W. Rose
Nicola Pozzobon
Ivan Reid
Date
2013, Jul 18

Definition at line 43 of file TTStubBuilder.h.

Constructor & Destructor Documentation

template<typename T >
TTStubBuilder< T >::TTStubBuilder ( const edm::ParameterSet iConfig)
explicit

Constructor.

Close class.

Implementation of methods

Here, in the header file, the methods which do not depend on the specific type <T> that can fit the template. Other methods, with type-specific features, are implemented in the source file.Constructors

Definition at line 109 of file TTStubBuilder.h.

References TTStubBuilder< T >::applyFE, TTStubBuilder< T >::clustersToken, TTStubBuilder< T >::ForbidMultipleStubs, edm::ParameterSet::getParameter(), TTStubBuilder< T >::high_rate_max_ring, TTStubBuilder< T >::maxStubs_2S, TTStubBuilder< T >::maxStubs_2S_CIC_5, TTStubBuilder< T >::maxStubs_PS, TTStubBuilder< T >::maxStubs_PS_CIC_10, TTStubBuilder< T >::maxStubs_PS_CIC_5, TTStubBuilder< T >::tedd1_maxring, and TTStubBuilder< T >::tedd2_maxring.

109  {
110  clustersToken = consumes<edmNew::DetSetVector<TTCluster<T> > >(iConfig.getParameter<edm::InputTag>("TTClusters"));
111  ForbidMultipleStubs = iConfig.getParameter<bool>("OnlyOnePerInputCluster");
112  applyFE = iConfig.getParameter<bool>("FEineffs");
113  maxStubs_2S = iConfig.getParameter<uint32_t>("CBClimit");
114  maxStubs_PS = iConfig.getParameter<uint32_t>("MPAlimit");
115  maxStubs_2S_CIC_5 = iConfig.getParameter<uint32_t>("SS5GCIClimit");
116  maxStubs_PS_CIC_5 = iConfig.getParameter<uint32_t>("PS5GCIClimit");
117  maxStubs_PS_CIC_10 = iConfig.getParameter<uint32_t>("PS10GCIClimit");
118  tedd1_maxring = iConfig.getParameter<uint32_t>("TEDD1Max10GRing");
119  tedd2_maxring = iConfig.getParameter<uint32_t>("TEDD2Max10GRing");
120  produces<edmNew::DetSetVector<TTCluster<T> > >("ClusterAccepted");
121  produces<edmNew::DetSetVector<TTStub<T> > >("StubAccepted");
122  produces<edmNew::DetSetVector<TTStub<T> > >("StubRejected");
123 
129 }
T getParameter(std::string const &) const
unsigned int maxStubs_PS_CIC_10
Definition: TTStubBuilder.h:77
unsigned int high_rate_max_ring[5]
Definition: TTStubBuilder.h:96
bool ForbidMultipleStubs
Definition: TTStubBuilder.h:55
unsigned int maxStubs_2S_CIC_5
Definition: TTStubBuilder.h:75
unsigned int tedd2_maxring
Definition: TTStubBuilder.h:80
unsigned int maxStubs_2S
Definition: TTStubBuilder.h:73
edm::EDGetTokenT< edmNew::DetSetVector< TTCluster< T > > > clustersToken
Definition: TTStubBuilder.h:54
unsigned int maxStubs_PS
Definition: TTStubBuilder.h:74
unsigned int maxStubs_PS_CIC_5
Definition: TTStubBuilder.h:76
unsigned int tedd1_maxring
Definition: TTStubBuilder.h:79
template<typename T >
TTStubBuilder< T >::~TTStubBuilder ( )
override

Destructor;.

Destructor.

Definition at line 133 of file TTStubBuilder.h.

133 {}

Member Function Documentation

template<typename T >
void TTStubBuilder< T >::beginRun ( const edm::Run run,
const edm::EventSetup iSetup 
)
overrideprivate

Mandatory methods.

Begin run.

Get the stub finding algorithm

Definition at line 137 of file TTStubBuilder.h.

References edm::EventSetup::get(), TTStubBuilder< T >::ievt, TTStubBuilder< T >::moduleStubs_CBC, TTStubBuilder< T >::moduleStubs_CIC, TTStubBuilder< T >::moduleStubs_MPA, and TTStubBuilder< T >::theStubFindingAlgoHandle.

137  {
140  ievt = 0;
141  moduleStubs_CIC.clear();
142  moduleStubs_MPA.clear();
143  moduleStubs_CBC.clear();
144 }
edm::ESHandle< TTStubAlgorithm< T > > theStubFindingAlgoHandle
Data members.
Definition: TTStubBuilder.h:53
Class to store the TTStubAlgorithm used in TTStubBuilder.
std::unordered_map< int, int > moduleStubs_MPA
Definition: TTStubBuilder.h:87
std::unordered_map< int, int > moduleStubs_CBC
Definition: TTStubBuilder.h:88
T get() const
Definition: EventSetup.h:73
std::unordered_map< int, std::vector< TTStub< Ref_Phase2TrackerDigi_ > > > moduleStubs_CIC
Temporary storage for stubs before max check.
Definition: TTStubBuilder.h:86
template<typename T >
void TTStubBuilder< T >::endRun ( const edm::Run run,
const edm::EventSetup iSetup 
)
overrideprivate

End run.

Definition at line 148 of file TTStubBuilder.h.

148 {}
template<>
void TTStubBuilder< Ref_Phase2TrackerDigi_ >::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
private

Implementation of methods of TTClusterBuilder.h.

Here, in the source file, the methods which do depend on the specific type <T> that can fit the template.

Author
Andrew W. Rose
Nicola Pozzobon
Ivan Reid
Date
2013, Jul 18Implement the producer

Prepare output

Get the Clusters already stored away

Go on only if both detectors have Clusters

Get the DetSets of the Clusters

If there are Clusters in both sensors you can try and make a Stub This is ~redundant

Create the vectors of objects to be passed to the FastFillers

Get chip size information

SV 21/11/17: tracker topology method should be updated, currently provide wrong nums

Loop over pairs of Clusters

Temporary storage to allow only one stub per inner cluster if requested in cfi

Build a temporary Stub

Check for compatibility

If the Stub is above threshold

Stub accepted

End of loop over outer clusters

Here tempOutput stores all the stubs from this inner cluster Check if there is need to store only one (if only one already, skip this step)

If so, sort the stubs by bend and keep only the first one (smallest bend)

Get to the second element (the switch above ensures there are min 2)

tempIter points now to the second element

Delete all-but-the first one from tempOutput

Here, tempOutput is either of size 1 (if entering the switch) either of size N with all the valid combinations ...

Now loop over the accepted stubs (1 or N) for this inner cluster

Get the stub

Put in the output

This means that ALL stubs go into the output

This means that only some of them do Put in the temporary output

Find out which MPA/CBC ASIC

Find out which CIC ASIC

Already a stub for this ASIC?

No, so new entry

Already a stub for this ASIC?

No, so new entry

Already a stub for this ASIC?

Already a stub for this ASIC?

No, so new entry

Already a stub for this ASIC?

No, so new entry

Already a stub for this ASIC?

Sort them by |bend| and pick up only the first N.

End of check on max number of stubs per module

End of nested loop

End of loop over pairs of Clusters

Create the FastFillers

End of loop over detector elements

Put output in the event (1) Get also the OrphanHandle of the accepted clusters

Now, correctly reset the output

Get the DetId and prepare the FastFiller

detid of the two components. This should be done via a TrackerTopology method that is not yet available.

Get the DetSets of the clusters

Get the DetSet of the stubs

Prepare the new DetSet to replace the current one Loop over the stubs

Create a temporary stub

Compare the clusters stored in the stub with the ones of this module

If no compatible clusters were found, skip to the next one

getter is in FULL-strip units, setter is in HALF-strip units

getter is in FULL-strip units, setter is in HALF-strip units

End of loop over stubs of this module

End of loop over stub DetSetVector

Put output in the event (2)

Definition at line 16 of file TTStubBuilder.cc.

References TTStub< T >::addClusterRef(), edmNew::DetSet< T >::begin(), edmNew::DetSetVector< T >::begin(), TTStub< T >::bendBE(), TTStub< T >::bendOffset(), TTStub< T >::clusterRef(), constexpr, TrackerGeometry::dets(), edmNew::DetSet< T >::empty(), edmNew::DetSetVector< T >::empty(), edmNew::DetSet< T >::end(), edmNew::DetSetVector< T >::end(), dqmdumpme::first, edm::EventSetup::get(), edm::Event::getByToken(), TrackerGeometry::getDetectorType(), TTStub< T >::getDetId(), mps_fire::i, edmNew::DetSetVector< T >::id(), TTStub< T >::innerClusterPosition(), createfilelist::int, TrackerTopology::isLower(), TrackerTopology::layer(), visualization-live-secondInstance_cfg::m, edmNew::makeRefTo(), TTStub< T >::moduleTypePS(), eostools::move(), TrackerTopology::partnerDetId(), TrackerGeometry::Ph2PSP, edm::ESHandle< T >::product(), edmNew::DetSetVector< T >::push_back(), edm::Event::put(), TTStub< T >::rawBend(), TTStub< T >::setBendBE(), TTStub< T >::setBendOffset(), TTStub< T >::setModuleTypePS(), TTStub< T >::setRawBend(), findQualityFiles::size, TrackerTopology::stack(), DetId::subdetId(), StripSubdetector::TID, TrackerTopology::tidRing(), TrackerTopology::tidWheel(), StripSubdetector::TOB, and PV2DBase< T, PVType, FrameType >::y().

16  {
17  //Retrieve tracker topology from geometry
19  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
20  const TrackerTopology* const tTopo = tTopoHandle.product();
22  iSetup.get<TrackerDigiGeometryRecord>().get(tGeomHandle);
23  const TrackerGeometry* const theTrackerGeom = tGeomHandle.product();
24 
26  auto ttClusterDSVForOutput = std::make_unique<edmNew::DetSetVector<TTCluster<Ref_Phase2TrackerDigi_>>>();
27  auto ttStubDSVForOutputTemp = std::make_unique<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>>();
28  auto ttStubDSVForOutputAccepted = std::make_unique<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>>();
29  auto ttStubDSVForOutputRejected = std::make_unique<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>>();
30 
31  static constexpr int CBCFailOffset = 500;
32  static constexpr int CICFailOffset = 1000;
33 
36  iEvent.getByToken(clustersToken, clusterHandle);
37 
38  int nmod = -1;
39 
40  // Loop over all the tracker elements
41 
42  // for (auto gd=theTrackerGeom->dets().begin(); gd != theTrackerGeom->dets().end(); gd++)
43  for (const auto& gd : theTrackerGeom->dets()) {
44  DetId detid = (*gd).geographicalId();
45  if (detid.subdetId() != StripSubdetector::TOB && detid.subdetId() != StripSubdetector::TID)
46  continue; // only run on OT
47  if (!tTopo->isLower(detid))
48  continue; // loop on the stacks: choose the lower arbitrarily
49  DetId lowerDetid = detid;
50  DetId upperDetid = tTopo->partnerDetId(detid);
51  DetId stackDetid = tTopo->stack(detid);
52  bool isPS = (theTrackerGeom->getDetectorType(stackDetid) == TrackerGeometry::ModuleType::Ph2PSP);
53 
54  bool is10G_PS = false;
55 
56  // Determine if this module is a 10G transmission scheme module
57  //
58  // sviret comment (221217): this info should be made available in conddb at some point
59  // not in TrackerTopology as some modules may switch between 10G and 5G transmission
60  // schemes during running period
61 
62  if (detid.subdetId() == StripSubdetector::TOB) {
63  if (tTopo->layer(detid) == 1)
64  is10G_PS = true;
65  } else if (detid.subdetId() == StripSubdetector::TID) {
66  if (tTopo->tidRing(detid) <= high_rate_max_ring[tTopo->tidWheel(detid) - 1])
67  is10G_PS = true;
68  }
69 
70  ++nmod;
71 
72  unsigned int maxStubs;
73  std::vector<std::pair<unsigned int, double>> bendMap;
74 
76  if (clusterHandle->find(lowerDetid) == clusterHandle->end() ||
77  clusterHandle->find(upperDetid) == clusterHandle->end())
78  continue;
79 
81  edmNew::DetSet<TTCluster<Ref_Phase2TrackerDigi_>> lowerClusters = (*clusterHandle)[lowerDetid];
82  edmNew::DetSet<TTCluster<Ref_Phase2TrackerDigi_>> upperClusters = (*clusterHandle)[upperDetid];
83 
87  if (lowerClusters.empty() || upperClusters.empty())
88  continue;
89 
91  std::vector<TTCluster<Ref_Phase2TrackerDigi_>> tempInner;
92  std::vector<TTCluster<Ref_Phase2TrackerDigi_>> tempOuter;
93  std::vector<TTStub<Ref_Phase2TrackerDigi_>> tempAccepted;
94  tempInner.clear();
95  tempOuter.clear();
96  tempAccepted.clear();
97 
99  int chipSize = 127;
100  if (isPS)
101  chipSize = 120;
102 
104  for (auto lowerClusterIter = lowerClusters.begin(); lowerClusterIter != lowerClusters.end(); ++lowerClusterIter) {
107  std::vector<TTStub<Ref_Phase2TrackerDigi_>> tempOutput;
108 
109  for (auto upperClusterIter = upperClusters.begin(); upperClusterIter != upperClusters.end(); ++upperClusterIter) {
111  TTStub<Ref_Phase2TrackerDigi_> tempTTStub(stackDetid);
112  tempTTStub.addClusterRef(edmNew::makeRefTo(clusterHandle, lowerClusterIter));
113  tempTTStub.addClusterRef(edmNew::makeRefTo(clusterHandle, upperClusterIter));
114  tempTTStub.setModuleTypePS(isPS);
115 
117  bool thisConfirmation = false;
118  int thisDisplacement = 999999;
119  int thisOffset = 0;
120  float thisHardBend = 0;
121 
122  theStubFindingAlgoHandle->PatternHitCorrelation(
123  thisConfirmation, thisDisplacement, thisOffset, thisHardBend, tempTTStub);
124  // Removed real offset. Ivan Reid 10/2019
125 
127  if (thisConfirmation) {
128  tempTTStub.setRawBend(thisDisplacement);
129  tempTTStub.setBendOffset(thisOffset);
130  tempTTStub.setBendBE(thisHardBend);
131  tempOutput.push_back(tempTTStub);
132  }
133  }
134 
137  if (ForbidMultipleStubs && tempOutput.size() > 1) {
139  std::sort(tempOutput.begin(), tempOutput.end(), TTStubBuilder<Ref_Phase2TrackerDigi_>::SortStubsBend);
140 
142  typename std::vector<TTStub<Ref_Phase2TrackerDigi_>>::iterator tempIter = tempOutput.begin();
143  ++tempIter;
144 
146 
148  tempOutput.erase(tempIter, tempOutput.end());
149  }
150 
153 
155  for (unsigned int iTempStub = 0; iTempStub < tempOutput.size(); ++iTempStub) {
157  const TTStub<Ref_Phase2TrackerDigi_>& tempTTStub = tempOutput[iTempStub];
158 
159  // A temporary stub, for FE problems
160  TTStub<Ref_Phase2TrackerDigi_> tempTTStub2(tempTTStub.getDetId());
161 
162  tempTTStub2.addClusterRef((tempTTStub.clusterRef(0)));
163  tempTTStub2.addClusterRef((tempTTStub.clusterRef(1)));
164  tempTTStub2.setRawBend(2. * tempTTStub.rawBend());
165  tempTTStub2.setBendOffset(2. * tempTTStub.bendOffset());
166  tempTTStub2.setBendBE(tempTTStub.bendBE());
167  tempTTStub2.setModuleTypePS(tempTTStub.moduleTypePS());
168 
170  if (!applyFE) // No dynamic inefficiencies (DEFAULT)
171  {
173  tempInner.push_back(*(tempTTStub.clusterRef(0)));
174  tempOuter.push_back(*(tempTTStub.clusterRef(1)));
175  tempAccepted.push_back(tempTTStub);
176  } else {
177  bool FEreject = false;
180  MeasurementPoint mp0 = tempTTStub.clusterRef(0)->findAverageLocalCoordinates();
181  int seg = static_cast<int>(mp0.y());
182  if (isPS)
183  seg = seg / 16;
184  int chip = 1000 * nmod + 10 * int(tempTTStub.innerClusterPosition() / chipSize) +
185  seg;
186  int CIC_chip = 10 * nmod + seg;
187 
188  // First look is the stub is passing trough the very front end (CBC/MPA)
189  (isPS) ? maxStubs = maxStubs_PS : maxStubs = maxStubs_2S;
190 
191  if (isPS) // MPA
192  {
193  if (moduleStubs_MPA.find(chip) == moduleStubs_MPA.end())
194  {
196  moduleStubs_MPA.emplace(chip, 1);
197  } else {
198  if (moduleStubs_MPA[chip] < int(maxStubs)) {
199  ++moduleStubs_MPA[chip];
200  } else {
201  FEreject = true;
202  }
203  }
204  } else // CBC
205  {
206  if (moduleStubs_CBC.find(chip) == moduleStubs_CBC.end())
207  {
209  moduleStubs_CBC.emplace(chip, 1);
210  } else {
211  if (moduleStubs_CBC[chip] < int(maxStubs)) {
212  ++moduleStubs_CBC[chip];
213  } else {
214  FEreject = true;
215  }
216  }
217  }
218 
219  // End of the MPA/CBC loop
220 
221  // If the stub has been already thrown out, there is no reason to include it into the CIC stream
222  // We keep is in the stub final container tough, but flagged as reject by FE
223 
224  if (FEreject) {
225  tempTTStub2.setRawBend(CBCFailOffset + 2. * tempTTStub.rawBend());
226  tempTTStub2.setBendOffset(CBCFailOffset + 2. * tempTTStub.bendOffset());
227 
228  tempInner.push_back(*(tempTTStub2.clusterRef(0)));
229  tempOuter.push_back(*(tempTTStub2.clusterRef(1)));
230  tempAccepted.push_back(tempTTStub2);
231  continue;
232  }
233 
234  maxStubs = isPS ? maxStubs_PS_CIC_5 : maxStubs_2S_CIC_5;
235 
236  if (is10G_PS)
237  maxStubs = maxStubs_PS_CIC_10;
238 
239  if (moduleStubs_CIC.find(CIC_chip) == moduleStubs_CIC.end())
240  {
241  if (moduleStubs_MPA.find(chip) == moduleStubs_MPA.end())
242  {
244  moduleStubs_MPA.emplace(chip, 1);
245  } else {
246  if (moduleStubs_MPA[chip] < int(maxStubs)) {
247  ++moduleStubs_MPA[chip];
248  } else {
249  FEreject = true;
250  }
251  }
252  } else // CBC
253  {
254  if (moduleStubs_CBC.find(chip) == moduleStubs_CBC.end())
255  {
257  moduleStubs_CBC.emplace(chip, 1);
258  } else {
259  if (moduleStubs_CBC[chip] < int(maxStubs)) {
260  ++moduleStubs_CBC[chip];
261  } else {
262  FEreject = true;
263  }
264  }
265  }
266 
267  // End of the MPA/CBC loop
268 
269  // If the stub has been already thrown out, there is no reason to include it into the CIC stream
270  // We keep is in the stub final container tough, but flagged as reject by FE
271 
272  if (FEreject) {
273  tempTTStub2.setRawBend(CBCFailOffset + 2. * tempTTStub.rawBend());
274  tempTTStub2.setBendOffset(CBCFailOffset + 2. * tempTTStub.bendOffset());
275 
276  tempInner.push_back(*(tempTTStub2.clusterRef(0)));
277  tempOuter.push_back(*(tempTTStub2.clusterRef(1)));
278  tempAccepted.push_back(tempTTStub2);
279  continue;
280  }
281 
282  (isPS) ? maxStubs = maxStubs_PS_CIC_5 : maxStubs = maxStubs_2S_CIC_5;
283 
284  if (is10G_PS)
285  maxStubs = maxStubs_PS_CIC_10;
286 
287  if (moduleStubs_CIC.find(CIC_chip) == moduleStubs_CIC.end())
288  {
289  bool CIC_reject = true;
290 
291  if (moduleStubs_CIC[CIC_chip].size() < maxStubs) {
292  moduleStubs_CIC[CIC_chip].push_back(tempTTStub); //We add the new stub
293  tempInner.push_back(*(tempTTStub.clusterRef(0)));
294  tempOuter.push_back(*(tempTTStub.clusterRef(1)));
295  tempAccepted.push_back(tempTTStub); // The stub is added
296  } else {
297  moduleStubs_CIC[CIC_chip].push_back(tempTTStub); //We add the new stub
298 
300  bendMap.clear();
301  bendMap.reserve(moduleStubs_CIC[CIC_chip].size());
302 
303  for (unsigned int i = 0; i < moduleStubs_CIC[CIC_chip].size(); ++i) {
304  bendMap.emplace_back(i, moduleStubs_CIC[CIC_chip].at(i).bendFE());
305  }
306 
307  std::sort(bendMap.begin(), bendMap.end(), TTStubBuilder<Ref_Phase2TrackerDigi_>::SortStubBendPairs);
308 
309  // bendMap contains link over all the stubs included in moduleStubs_CIC[CIC_chip]
310 
311  for (unsigned int i = 0; i < maxStubs; ++i) {
312  // The stub we have added is among the first ones, add it
313  if (bendMap[i].first == moduleStubs_CIC[CIC_chip].size() - 1) {
314  CIC_reject = false;
315  }
316  }
317 
318  if (CIC_reject) // The stub added does not pass the cut
319  {
320  tempTTStub2.setRawBend(CICFailOffset + 2. * tempTTStub.rawBend());
321  tempTTStub2.setBendOffset(CICFailOffset + 2. * tempTTStub.bendOffset());
322 
323  tempInner.push_back(*(tempTTStub2.clusterRef(0)));
324  tempOuter.push_back(*(tempTTStub2.clusterRef(1)));
325  tempAccepted.push_back(tempTTStub2);
326  } else {
327  tempInner.push_back(*(tempTTStub.clusterRef(0)));
328  tempOuter.push_back(*(tempTTStub.clusterRef(1)));
329  tempAccepted.push_back(tempTTStub); // The stub is added
330  }
331  }
332  }
333  }
334  }
335  }
336 
338  if (!tempInner.empty()) {
339  typename edmNew::DetSetVector<TTCluster<Ref_Phase2TrackerDigi_>>::FastFiller lowerOutputFiller(
340  *ttClusterDSVForOutput, lowerDetid);
341  for (unsigned int m = 0; m < tempInner.size(); m++) {
342  lowerOutputFiller.push_back(tempInner.at(m));
343  }
344  if (lowerOutputFiller.empty())
345  lowerOutputFiller.abort();
346  }
347 
348  if (!tempOuter.empty()) {
349  typename edmNew::DetSetVector<TTCluster<Ref_Phase2TrackerDigi_>>::FastFiller upperOutputFiller(
350  *ttClusterDSVForOutput, upperDetid);
351  for (unsigned int m = 0; m < tempOuter.size(); m++) {
352  upperOutputFiller.push_back(tempOuter.at(m));
353  }
354  if (upperOutputFiller.empty())
355  upperOutputFiller.abort();
356  }
357 
358  if (!tempAccepted.empty()) {
359  typename edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>::FastFiller tempAcceptedFiller(
360  *ttStubDSVForOutputTemp, stackDetid);
361  for (unsigned int m = 0; m < tempAccepted.size(); m++) {
362  tempAcceptedFiller.push_back(tempAccepted.at(m));
363  }
364  if (tempAcceptedFiller.empty())
365  tempAcceptedFiller.abort();
366  }
367 
368  }
369 
373  iEvent.put(std::move(ttClusterDSVForOutput), "ClusterAccepted");
374 
376  typename edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>::const_iterator stubDetIter;
377 
378  for (stubDetIter = ttStubDSVForOutputTemp->begin(); stubDetIter != ttStubDSVForOutputTemp->end(); ++stubDetIter) {
380  DetId thisStackedDetId = stubDetIter->id();
381  typename edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>::FastFiller acceptedOutputFiller(
382  *ttStubDSVForOutputAccepted, thisStackedDetId);
383 
386  DetId lowerDetid = thisStackedDetId + 1;
387  DetId upperDetid = thisStackedDetId + 2;
388 
390  edmNew::DetSet<TTCluster<Ref_Phase2TrackerDigi_>> lowerClusters = (*ttClusterAcceptedHandle)[lowerDetid];
391  edmNew::DetSet<TTCluster<Ref_Phase2TrackerDigi_>> upperClusters = (*ttClusterAcceptedHandle)[upperDetid];
392 
394  edmNew::DetSet<TTStub<Ref_Phase2TrackerDigi_>> theseStubs = (*ttStubDSVForOutputTemp)[thisStackedDetId];
395 
398  typename edmNew::DetSet<TTCluster<Ref_Phase2TrackerDigi_>>::iterator clusterIter;
399  typename edmNew::DetSet<TTStub<Ref_Phase2TrackerDigi_>>::iterator stubIter;
400  for (stubIter = theseStubs.begin(); stubIter != theseStubs.end(); ++stubIter) {
402  TTStub<Ref_Phase2TrackerDigi_> tempTTStub(stubIter->getDetId());
403 
406  lowerClusterToBeReplaced = stubIter->clusterRef(0);
407  const edm::Ref<edmNew::DetSetVector<TTCluster<Ref_Phase2TrackerDigi_>>, TTCluster<Ref_Phase2TrackerDigi_>>&
408  upperClusterToBeReplaced = stubIter->clusterRef(1);
409 
410  bool lowerOK = false;
411  bool upperOK = false;
412 
413  for (clusterIter = lowerClusters.begin(); clusterIter != lowerClusters.end() && !lowerOK; ++clusterIter) {
414  if (clusterIter->getHits() == lowerClusterToBeReplaced->getHits()) {
415  tempTTStub.addClusterRef(edmNew::makeRefTo(ttClusterAcceptedHandle, clusterIter));
416  lowerOK = true;
417  }
418  }
419 
420  for (clusterIter = upperClusters.begin(); clusterIter != upperClusters.end() && !upperOK; ++clusterIter) {
421  if (clusterIter->getHits() == upperClusterToBeReplaced->getHits()) {
422  tempTTStub.addClusterRef(edmNew::makeRefTo(ttClusterAcceptedHandle, clusterIter));
423  upperOK = true;
424  }
425  }
426 
428  if (!lowerOK || !upperOK)
429  continue;
430 
431  tempTTStub.setRawBend(2. * stubIter->rawBend());
432  tempTTStub.setBendOffset(
433  2. * stubIter->bendOffset());
434  tempTTStub.setBendBE(stubIter->bendBE());
435  tempTTStub.setModuleTypePS(stubIter->moduleTypePS());
436 
437  acceptedOutputFiller.push_back(tempTTStub);
438 
439  }
440 
441  if (acceptedOutputFiller.empty())
442  acceptedOutputFiller.abort();
443 
444  }
445 
447  iEvent.put(std::move(ttStubDSVForOutputAccepted), "StubAccepted");
448  iEvent.put(std::move(ttStubDSVForOutputRejected), "StubRejected");
449 
450  ++ievt;
451  if (ievt % 8 == 0)
452  moduleStubs_CIC.clear(); // Everything is cleared up after 8BX
453  if (ievt % 2 == 0)
454  moduleStubs_MPA.clear(); // Everything is cleared up after 2BX
455  moduleStubs_CBC.clear(); // Everything is cleared up after everyBX
456 }
size
Write out results.
edm::Ref< typename HandleT::element_type, typename HandleT::element_type::value_type::value_type > makeRefTo(const HandleT &iHandle, typename HandleT::element_type::value_type::const_iterator itIter)
const_iterator end(bool update=false) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
unsigned int maxStubs_PS_CIC_10
Definition: TTStubBuilder.h:77
T y() const
Definition: PV2DBase.h:44
void setBendOffset(int anOffset)
Definition: TTStub.h:180
unsigned int tidRing(const DetId &id) const
edm::ESHandle< TTStubAlgorithm< T > > theStubFindingAlgoHandle
Data members.
Definition: TTStubBuilder.h:53
void setModuleTypePS(bool isPSModule)
set whether this is a PS module or not;
Definition: TTStub.h:190
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
double rawBend() const
Trigger info.
Definition: TTStub.h:165
void setBendBE(float aBend)
setBendBE(): In HALF-STRIP units! Reduced resolution in BE boards. Rename of setHardwareBend() ...
Definition: TTStub.h:185
unsigned int high_rate_max_ring[5]
Definition: TTStubBuilder.h:96
double innerClusterPosition() const
Definition: TTStub.h:198
unsigned int tidWheel(const DetId &id) const
DetId getDetId() const
Detector element.
Definition: TTStub.h:44
double bendOffset() const
Definition: TTStub.h:175
bool ForbidMultipleStubs
Definition: TTStubBuilder.h:55
id_type id(size_t cell) const
DetId partnerDetId(const DetId &id) const
bool isLower(const DetId &id) const
unsigned int maxStubs_2S_CIC_5
Definition: TTStubBuilder.h:75
const DetContainer & dets() const override
Returm a vector of all GeomDet (including all GeomDetUnits)
void setRawBend(int aDisplacement)
Definition: TTStub.h:170
bool empty() const
Definition: DetSetNew.h:70
unsigned int maxStubs_2S
Definition: TTStubBuilder.h:73
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
edm::EDGetTokenT< edmNew::DetSetVector< TTCluster< T > > > clustersToken
Definition: TTStubBuilder.h:54
static constexpr auto TOB
ModuleType getDetectorType(DetId) const
Class to store the L1 Track Trigger stubs.
Definition: TTStub.h:22
std::unordered_map< int, int > moduleStubs_MPA
Definition: TTStubBuilder.h:87
unsigned int maxStubs_PS
Definition: TTStubBuilder.h:74
std::unordered_map< int, int > moduleStubs_CBC
Definition: TTStubBuilder.h:88
const edm::Ref< edmNew::DetSetVector< TTCluster< T > >, TTCluster< T > > & clusterRef(unsigned int hitStackMember) const
Clusters composing the Stub – see https://twiki.cern.ch/twiki/bin/viewauth/CMS/SLHCTrackerTriggerSWT...
Definition: TTStub.h:150
Definition: DetId.h:17
NOTE: this is needed even if it seems not.
Definition: TTCluster.h:27
uint32_t stack(const DetId &id) const
void addClusterRef(edm::Ref< edmNew::DetSetVector< TTCluster< T > >, TTCluster< T > > aTTCluster)
Add a cluster reference, depending on which stack member it is on (inner = 0, outer = 1) ...
Definition: TTStub.h:156
unsigned int layer(const DetId &id) const
double bendBE() const
BendBE(): In FULL-STRIP units! Reduced resolution from BE boards. Rename of getHardwareBend() ...
Definition: TTStub.h:211
unsigned int maxStubs_PS_CIC_5
Definition: TTStubBuilder.h:76
iterator end()
Definition: DetSetNew.h:56
T get() const
Definition: EventSetup.h:73
std::unordered_map< int, std::vector< TTStub< Ref_Phase2TrackerDigi_ > > > moduleStubs_CIC
Temporary storage for stubs before max check.
Definition: TTStubBuilder.h:86
static constexpr auto TID
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:511
const_iterator begin(bool update=false) const
#define constexpr
Plugin to load the Stub finding algorithm and produce the collection of Stubs that goes in the event ...
Definition: TTStubBuilder.h:43
bool moduleTypePS() const
check if a PS module
Definition: TTStub.h:194
iterator begin()
Definition: DetSetNew.h:54
template<typename T >
void TTStubBuilder< T >::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate
template<>
void TTStubBuilder< Ref_Phase2TrackerDigi_ >::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
private

Implement the producer.

template<typename T >
bool TTStubBuilder< T >::SortStubBendPairs ( const std::pair< unsigned int, double > &  left,
const std::pair< unsigned int, double > &  right 
)
staticprivate

Sort routine for stub ordering.

Sorting method for stubs NOTE: this must be static!

Definition at line 152 of file TTStubBuilder.h.

References funct::abs().

153  {
154  return std::abs(left.second) < std::abs(right.second);
155 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
template<typename T >
bool TTStubBuilder< T >::SortStubsBend ( const TTStub< T > &  left,
const TTStub< T > &  right 
)
staticprivate

Analogous sorting routine directly from stubs.

Definition at line 159 of file TTStubBuilder.h.

References funct::abs(), TTStub< T >::bendFE(), iEvent, and TTStubBuilder< T >::produce().

159  {
160  return std::abs(left.bendFE()) < std::abs(right.bendFE());
161 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double bendFE() const
BendFE(): In FULL-STRIP units from FE! Rename of getTriggerBend()
Definition: TTStub.h:203

Member Data Documentation

template<typename T >
bool TTStubBuilder< T >::applyFE
private

Definition at line 70 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::TTStubBuilder().

template<typename T >
edm::EDGetTokenT<edmNew::DetSetVector<TTCluster<T> > > TTStubBuilder< T >::clustersToken
private

Definition at line 54 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::TTStubBuilder().

template<typename T >
bool TTStubBuilder< T >::ForbidMultipleStubs
private

Definition at line 55 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::TTStubBuilder().

template<typename T >
unsigned int TTStubBuilder< T >::high_rate_max_ring[5]
private

Definition at line 96 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::TTStubBuilder().

template<typename T >
int TTStubBuilder< T >::ievt
private

Definition at line 82 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::beginRun().

template<typename T >
unsigned int TTStubBuilder< T >::maxStubs_2S
private

Definition at line 73 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::TTStubBuilder().

template<typename T >
unsigned int TTStubBuilder< T >::maxStubs_2S_CIC_5
private

Definition at line 75 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::TTStubBuilder().

template<typename T >
unsigned int TTStubBuilder< T >::maxStubs_PS
private

Definition at line 74 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::TTStubBuilder().

template<typename T >
unsigned int TTStubBuilder< T >::maxStubs_PS_CIC_10
private

Definition at line 77 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::TTStubBuilder().

template<typename T >
unsigned int TTStubBuilder< T >::maxStubs_PS_CIC_5
private

Definition at line 76 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::TTStubBuilder().

template<typename T >
std::unordered_map<int, int> TTStubBuilder< T >::moduleStubs_CBC
private

Definition at line 88 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::beginRun().

template<typename T >
std::unordered_map<int, std::vector<TTStub<Ref_Phase2TrackerDigi_> > > TTStubBuilder< T >::moduleStubs_CIC
private

Temporary storage for stubs before max check.

Definition at line 86 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::beginRun().

template<typename T >
std::unordered_map<int, int> TTStubBuilder< T >::moduleStubs_MPA
private

Definition at line 87 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::beginRun().

template<typename T >
unsigned int TTStubBuilder< T >::tedd1_maxring
private

Definition at line 79 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::TTStubBuilder().

template<typename T >
unsigned int TTStubBuilder< T >::tedd2_maxring
private

Definition at line 80 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::TTStubBuilder().

template<typename T >
edm::ESHandle<TTStubAlgorithm<T> > TTStubBuilder< T >::theStubFindingAlgoHandle
private

Data members.

Definition at line 53 of file TTStubBuilder.h.

Referenced by TTStubBuilder< T >::beginRun().