Main Page
Namespaces
Classes
Package Documentation
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Pages
src
RecoTracker
TkSeedingLayers
src
HitExtractorSTRP.h
Go to the documentation of this file.
1
#ifndef RecoTracker_TkSeedingLayers_HitExtractorSTRP_H
2
#define RecoTracker_TkSeedingLayers_HitExtractorSTRP_H
3
4
#include "
RecoTracker/TkSeedingLayers/interface/SeedingLayer.h
"
5
#include "
FWCore/Utilities/interface/InputTag.h
"
6
#include "
HitExtractor.h
"
7
8
#include "
DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h
"
9
10
#include <vector>
11
class
DetLayer
;
12
13
namespace
edm {
14
template
<
typename
T>
class
ContainerMask
;
15
}
16
17
namespace
ctfseeding {
18
19
class
HitExtractorSTRP
:
public
HitExtractor
{
20
21
public
:
22
typedef
SiStripRecHit2D::ClusterRef
SiStripClusterRef
;
23
24
HitExtractorSTRP
(
const
DetLayer
* detLayer,
SeedingLayer::Side
& side,
int
idLayer);
25
virtual
~HitExtractorSTRP
(){}
26
27
virtual
HitExtractor::Hits
hits
(
const
SeedingLayer
& sl,
const
edm::Event
& ,
const
edm::EventSetup
& )
const
;
28
virtual
HitExtractorSTRP
*
clone
()
const
{
return
new
HitExtractorSTRP
(*
this
); }
29
30
void
useMatchedHits
(
const
edm::InputTag
&
m
) {
hasMatchedHits
=
true
;
theMatchedHits
=
m
; }
31
void
useRPhiHits
(
const
edm::InputTag
&
m
) {
hasRPhiHits
=
true
;
theRPhiHits
=
m
; }
32
void
useStereoHits
(
const
edm::InputTag
&
m
) {
hasStereoHits
=
true
;
theStereoHits
=
m
; }
33
void
useRingSelector
(
int
minRing,
int
maxRing);
34
void
useSimpleRphiHitsCleaner
(
bool
use) {
hasSimpleRphiHitsCleaner
= use;}
35
36
void
cleanedOfClusters
(
const
edm::Event
& ev,
HitExtractor::Hits
&
hits
,
bool
matched,
unsigned
int
cleanFrom=0)
const
;
37
38
bool
skipThis
(
TransientTrackingRecHit::ConstRecHitPointer
& ptr,
edm::Handle
<
edm::ContainerMask
<
edmNew::DetSetVector<SiStripCluster>
> > & stripClusterMask,
39
TransientTrackingRecHit::ConstRecHitPointer
& replaceMe)
const
;
40
41
bool
skipThis
(
OmniClusterRef
const
& clus,
edm::Handle
<
edm::ContainerMask
<
edmNew::DetSetVector<SiStripCluster>
> > & stripClusterMask)
const
;
42
43
void
project
(
TransientTrackingRecHit::ConstRecHitPointer
& ptr,
44
const
SiStripRecHit2D *
hit
,
45
TransientTrackingRecHit::ConstRecHitPointer
& replaceMe)
const
;
46
void
setNoProjection
()
const
{
failProjection
=
true
;};
47
private
:
48
bool
ringRange
(
int
ring
)
const
;
49
private
:
50
const
DetLayer
*
theLayer
;
51
SeedingLayer::Side
theSide
;
52
mutable
const
SeedingLayer
*
theSLayer
;
53
int
theIdLayer
;
54
bool
hasMatchedHits
;
edm::InputTag
theMatchedHits
;
55
bool
hasRPhiHits
;
edm::InputTag
theRPhiHits
;
56
bool
hasStereoHits
;
edm::InputTag
theStereoHits
;
57
bool
hasRingSelector
;
int
theMinRing
,
theMaxRing
;
58
bool
hasSimpleRphiHitsCleaner
;
59
mutable
bool
failProjection
;
60
};
61
62
}
63
#endif
ctfseeding::HitExtractorSTRP::hasSimpleRphiHitsCleaner
bool hasSimpleRphiHitsCleaner
Definition:
HitExtractorSTRP.h:58
ctfseeding::HitExtractorSTRP::HitExtractorSTRP
HitExtractorSTRP(const DetLayer *detLayer, SeedingLayer::Side &side, int idLayer)
Definition:
HitExtractorSTRP.cc:27
ctfseeding::HitExtractorSTRP::useMatchedHits
void useMatchedHits(const edm::InputTag &m)
Definition:
HitExtractorSTRP.h:30
edm::ContainerMask
Definition:
ContainerMask.h:36
ConstReferenceCountingPointer< TransientTrackingRecHit >
ctfseeding::HitExtractorSTRP::theMaxRing
int theMaxRing
Definition:
HitExtractorSTRP.h:57
ctfseeding::HitExtractorSTRP::ringRange
bool ringRange(int ring) const
Definition:
HitExtractorSTRP.cc:41
HitExtractor.h
ctfseeding::HitExtractor
Definition:
HitExtractor.h:14
ctfseeding::HitExtractorSTRP::theSide
SeedingLayer::Side theSide
Definition:
HitExtractorSTRP.h:51
edm::Handle
Definition:
AssociativeIterator.h:48
SiStripRecHit2D.h
ctfseeding::SeedingLayer
Definition:
SeedingLayer.h:18
ctfseeding::HitExtractorSTRP::useRPhiHits
void useRPhiHits(const edm::InputTag &m)
Definition:
HitExtractorSTRP.h:31
ContainerMask
ctfseeding::HitExtractorSTRP::theStereoHits
edm::InputTag theStereoHits
Definition:
HitExtractorSTRP.h:56
ctfseeding::HitExtractorSTRP::hasRPhiHits
bool hasRPhiHits
Definition:
HitExtractorSTRP.h:55
ctfseeding::HitExtractorSTRP::hits
virtual HitExtractor::Hits hits(const SeedingLayer &sl, const edm::Event &, const edm::EventSetup &) const
Definition:
HitExtractorSTRP.cc:127
ctfseeding::HitExtractorSTRP::hasStereoHits
bool hasStereoHits
Definition:
HitExtractorSTRP.h:56
ctfseeding::HitExtractorSTRP::theMinRing
int theMinRing
Definition:
HitExtractorSTRP.h:57
edmNew::DetSetVector
Definition:
DetSetNew.h:8
edm::EventSetup
Definition:
EventSetup.h:44
DetLayer
Definition:
DetLayer.h:26
ctfseeding::HitExtractorSTRP::hasMatchedHits
bool hasMatchedHits
Definition:
HitExtractorSTRP.h:54
ctfseeding::HitExtractorSTRP::setNoProjection
void setNoProjection() const
Definition:
HitExtractorSTRP.h:46
relativeConstraints.ring
tuple ring
Definition:
relativeConstraints.py:67
ctfseeding::HitExtractorSTRP::theSLayer
const SeedingLayer * theSLayer
Definition:
HitExtractorSTRP.h:52
ctfseeding::HitExtractorSTRP
Definition:
HitExtractorSTRP.h:19
ctfseeding::HitExtractorSTRP::useRingSelector
void useRingSelector(int minRing, int maxRing)
Definition:
HitExtractorSTRP.cc:34
m
int m
Definition:
DTDataIntegrityTask.cc:33
ctfseeding::HitExtractorSTRP::cleanedOfClusters
void cleanedOfClusters(const edm::Event &ev, HitExtractor::Hits &hits, bool matched, unsigned int cleanFrom=0) const
Definition:
HitExtractorSTRP.cc:96
ctfseeding::HitExtractorSTRP::theLayer
const DetLayer * theLayer
Definition:
HitExtractorSTRP.h:50
ctfseeding::HitExtractorSTRP::theRPhiHits
edm::InputTag theRPhiHits
Definition:
HitExtractorSTRP.h:55
ctfseeding::HitExtractorSTRP::clone
virtual HitExtractorSTRP * clone() const
Definition:
HitExtractorSTRP.h:28
ctfseeding::HitExtractorSTRP::useSimpleRphiHitsCleaner
void useSimpleRphiHitsCleaner(bool use)
Definition:
HitExtractorSTRP.h:34
ctfseeding::HitExtractorSTRP::project
void project(TransientTrackingRecHit::ConstRecHitPointer &ptr, const SiStripRecHit2D *hit, TransientTrackingRecHit::ConstRecHitPointer &replaceMe) const
Definition:
HitExtractorSTRP.cc:56
OmniClusterRef
Definition:
OmniClusterRef.h:12
ctfseeding::HitExtractorSTRP::skipThis
bool skipThis(TransientTrackingRecHit::ConstRecHitPointer &ptr, edm::Handle< edm::ContainerMask< edmNew::DetSetVector< SiStripCluster > > > &stripClusterMask, TransientTrackingRecHit::ConstRecHitPointer &replaceMe) const
Definition:
HitExtractorSTRP.cc:67
SeedingLayer.h
hit
Definition:
SiStripHitEffFromCalibTree.cc:87
edm::InputTag
Definition:
InputTag.h:12
InputTag.h
ctfseeding::HitExtractorSTRP::theIdLayer
int theIdLayer
Definition:
HitExtractorSTRP.h:53
ctfseeding::HitExtractorSTRP::SiStripClusterRef
SiStripRecHit2D::ClusterRef SiStripClusterRef
Definition:
HitExtractorSTRP.h:22
ctfseeding::HitExtractorSTRP::~HitExtractorSTRP
virtual ~HitExtractorSTRP()
Definition:
HitExtractorSTRP.h:25
ctfseeding::HitExtractorSTRP::useStereoHits
void useStereoHits(const edm::InputTag &m)
Definition:
HitExtractorSTRP.h:32
edm::Event
Definition:
Event.h:50
ctfseeding::SeedingLayer::Side
Side
Definition:
SeedingLayer.h:20
ctfseeding::HitExtractorSTRP::failProjection
bool failProjection
Definition:
HitExtractorSTRP.h:59
ctfseeding::HitExtractorSTRP::hasRingSelector
bool hasRingSelector
Definition:
HitExtractorSTRP.h:57
ctfseeding::HitExtractorSTRP::theMatchedHits
edm::InputTag theMatchedHits
Definition:
HitExtractorSTRP.h:54
ctfseeding::HitExtractor::Hits
std::vector< TransientTrackingRecHit::ConstRecHitPointer > Hits
Definition:
HitExtractor.h:16
Generated for CMSSW Reference Manual by
1.8.5