Main Page
Namespaces
Classes
Package Documentation
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Pages
RecoLocalTracker
SiPixelRecHits
interface
SiPixelRecHitConverter.h
Go to the documentation of this file.
1
#ifndef RecoLocalTracker_SiPixelRecHits_SiPixelRecHitConverter_h
2
#define RecoLocalTracker_SiPixelRecHits_SiPixelRecHitConverter_h
3
4
//---------------------------------------------------------------------------
32
//---------------------------------------------------------------------------
33
34
//--- Base class for CPEs:
35
#include "
RecoLocalTracker/ClusterParameterEstimator/interface/PixelClusterParameterEstimator.h
"
36
//&&& #include "RecoLocalTracker/SiPixelRecHits/interface/PixelCPEBase.h"
37
38
//--- Geometry + DataFormats
39
#include "
Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h
"
40
#include "
DataFormats/SiPixelCluster/interface/SiPixelCluster.h
"
41
#include "
DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h
"
42
#include "
DataFormats/Common/interface/DetSetVector.h
"
43
44
//--- Framework
45
#include "
FWCore/Framework/interface/EDProducer.h
"
46
#include "
FWCore/Framework/interface/Event.h
"
47
#include "
FWCore/Framework/interface/EventSetup.h
"
48
49
//#define TP_OLD
50
#ifdef TP_OLD
51
#include "FWCore/Framework/interface/Handle.h"
52
#else
53
#include "
DataFormats/Common/interface/Handle.h
"
54
#endif
55
#include "
FWCore/Framework/interface/ESHandle.h
"
56
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
57
#include "
FWCore/Utilities/interface/InputTag.h
"
58
59
class
MagneticField
;
60
namespace
cms
61
{
62
class
SiPixelRecHitConverter
:
public
edm::EDProducer
63
{
64
public
:
65
//--- Constructor, virtual destructor (just in case)
66
explicit
SiPixelRecHitConverter
(
const
edm::ParameterSet
&
conf
);
67
virtual
~SiPixelRecHitConverter
();
68
69
//--- Factory method to make CPE's depending on the ParameterSet
70
//--- Not sure if we need to make more than one CPE to run concurrently
71
//--- on different parts of the detector (e.g., one for the barrel and the
72
//--- one for the forward). The way the CPE's are written now, it's
73
//--- likely we can use one (and they will switch internally), or
74
//--- make two of the same but configure them differently. We need a more
75
//--- realistic use case...
76
77
//--- The top-level event method.
78
virtual
void
produce
(
edm::Event
&
e
,
const
edm::EventSetup
&
c
);
79
80
// Begin Job
81
//virtual void beginJob();
82
virtual
void
beginJob
();
83
84
//--- Execute the position estimator algorithm(s).
85
//--- New interface with DetSetVector
86
void
run
(
const
edmNew::DetSetVector<SiPixelCluster>
&
input
,
87
SiPixelRecHitCollectionNew
&
output
,
88
edm::ESHandle<TrackerGeometry>
&
geom
);
89
90
void
run
(
edm::Handle
<
edmNew::DetSetVector<SiPixelCluster>
> inputhandle,
91
SiPixelRecHitCollectionNew
& output,
92
edm::ESHandle<TrackerGeometry>
& geom);
93
94
private
:
95
edm::ParameterSet
conf_
;
96
// TO DO: maybe allow a map of pointers?
97
std::string
cpeName_
;
// what the user said s/he wanted
98
const
PixelClusterParameterEstimator
*
cpe_
;
// what we got (for now, one ptr to base class)
99
//&&& PixelCPEBase * cpe_; // what we got (for now, one ptr to base class)
100
bool
ready_
;
// needed CPE's valid => good to go!
101
edm::InputTag
src_
;
102
int
theVerboseLevel
;
// algorithm's verbosity
103
bool
m_newCont
;
// save also in emdNew::DetSetVector
104
};
105
}
106
107
108
#endif
cms::SiPixelRecHitConverter::src_
edm::InputTag src_
Definition:
SiPixelRecHitConverter.h:101
cms::SiPixelRecHitConverter
Definition:
SiPixelRecHitConverter.h:62
cms::SiPixelRecHitConverter::run
void run(const edmNew::DetSetVector< SiPixelCluster > &input, SiPixelRecHitCollectionNew &output, edm::ESHandle< TrackerGeometry > &geom)
Event.h
LaserDQM_cfg.input
tuple input
Definition:
LaserDQM_cfg.py:38
EventSetup.h
MagneticField
Definition:
MagneticField.h:18
cms::SiPixelRecHitConverter::cpeName_
std::string cpeName_
Definition:
SiPixelRecHitConverter.h:97
Handle.h
edm::Handle
Definition:
AssociativeIterator.h:48
cms::SiPixelRecHitConverter::ready_
bool ready_
Definition:
SiPixelRecHitConverter.h:100
edm::EDProducer
Definition:
EDProducer.h:21
DetSetVector.h
ParameterSet.h
cms::SiPixelRecHitConverter::SiPixelRecHitConverter
SiPixelRecHitConverter(const edm::ParameterSet &conf)
Constructor: set the ParameterSet and defer all thinking to setupCPE().
Definition:
SiPixelRecHitConverter.cc:48
PixelClusterParameterEstimator.h
edm::ESHandle< TrackerGeometry >
relativeConstraints.geom
list geom
Definition:
relativeConstraints.py:71
SiPixelRecHitCollection.h
ESHandle.h
edmNew::DetSetVector
Definition:
DetSetNew.h:8
edm::EventSetup
Definition:
EventSetup.h:44
cms::SiPixelRecHitConverter::beginJob
virtual void beginJob()
Definition:
SiPixelRecHitConverter.cc:71
cms::SiPixelRecHitConverter::conf_
edm::ParameterSet conf_
Definition:
SiPixelRecHitConverter.h:95
dbtoconf.conf
tuple conf
Definition:
dbtoconf.py:185
EDProducer.h
cms::SiPixelRecHitConverter::~SiPixelRecHitConverter
virtual ~SiPixelRecHitConverter()
Definition:
SiPixelRecHitConverter.cc:63
convertSQLitetoXML_cfg.output
tuple output
Definition:
convertSQLitetoXML_cfg.py:32
cms::SiPixelRecHitConverter::produce
virtual void produce(edm::Event &e, const edm::EventSetup &c)
The "Event" entrypoint: gets called by framework for every event.
Definition:
SiPixelRecHitConverter.cc:78
trackerHits.c
tuple c
Definition:
trackerHits.py:26
PixelClusterParameterEstimator
Definition:
PixelClusterParameterEstimator.h:9
alignCSCRings.e
list e
Definition:
alignCSCRings.py:90
edm::InputTag
Definition:
InputTag.h:12
InputTag.h
TrackerGeometry.h
cms::SiPixelRecHitConverter::m_newCont
bool m_newCont
Definition:
SiPixelRecHitConverter.h:103
edm::ParameterSet
Definition:
ParameterSet.h:35
edm::Event
Definition:
Event.h:50
SiPixelCluster.h
cms::SiPixelRecHitConverter::theVerboseLevel
int theVerboseLevel
Definition:
SiPixelRecHitConverter.h:102
cms::SiPixelRecHitConverter::cpe_
const PixelClusterParameterEstimator * cpe_
Definition:
SiPixelRecHitConverter.h:98
Generated for CMSSW Reference Manual by
1.8.5