CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFClusterProducer.cc
Go to the documentation of this file.
2 
3 #include <memory>
4 
6 
9 
13 
14 using namespace std;
15 using namespace edm;
16 
18 {
19 
20  verbose_ =
21  iConfig.getUntrackedParameter<bool>("verbose",false);
22 
23 
24 
25  // parameters for clustering
26 
27  double threshBarrel =
28  iConfig.getParameter<double>("thresh_Barrel");
29  double threshSeedBarrel =
30  iConfig.getParameter<double>("thresh_Seed_Barrel");
31 
32  double threshPtBarrel =
33  iConfig.getParameter<double>("thresh_Pt_Barrel");
34  double threshPtSeedBarrel =
35  iConfig.getParameter<double>("thresh_Pt_Seed_Barrel");
36 
37  double threshCleanBarrel =
38  iConfig.getParameter<double>("thresh_Clean_Barrel");
39  std::vector<double> minS4S1CleanBarrel =
40  iConfig.getParameter< std::vector<double> >("minS4S1_Clean_Barrel");
41 
42  double threshEndcap =
43  iConfig.getParameter<double>("thresh_Endcap");
44  double threshSeedEndcap =
45  iConfig.getParameter<double>("thresh_Seed_Endcap");
46 
47  double threshPtEndcap =
48  iConfig.getParameter<double>("thresh_Pt_Endcap");
49  double threshPtSeedEndcap =
50  iConfig.getParameter<double>("thresh_Pt_Seed_Endcap");
51 
52  double threshCleanEndcap =
53  iConfig.getParameter<double>("thresh_Clean_Endcap");
54  std::vector<double> minS4S1CleanEndcap =
55  iConfig.getParameter< std::vector<double> >("minS4S1_Clean_Endcap");
56 
57  double threshDoubleSpikeBarrel =
58  iConfig.getParameter<double>("thresh_DoubleSpike_Barrel");
59  double minS6S2DoubleSpikeBarrel =
60  iConfig.getParameter<double>("minS6S2_DoubleSpike_Barrel");
61  double threshDoubleSpikeEndcap =
62  iConfig.getParameter<double>("thresh_DoubleSpike_Endcap");
63  double minS6S2DoubleSpikeEndcap =
64  iConfig.getParameter<double>("minS6S2_DoubleSpike_Endcap");
65 
66  int nNeighbours =
67  iConfig.getParameter<int>("nNeighbours");
68 
69 // double posCalcP1 =
70 // iConfig.getParameter<double>("posCalcP1");
71 
72  int posCalcNCrystal =
73  iConfig.getParameter<int>("posCalcNCrystal");
74 
75  double showerSigma =
76  iConfig.getParameter<double>("showerSigma");
77 
78  bool useCornerCells =
79  iConfig.getParameter<bool>("useCornerCells");
80 
81  bool cleanRBXandHPDs =
82  iConfig.getParameter<bool>("cleanRBXandHPDs");
83 
84 
85  clusterAlgo_.setThreshBarrel( threshBarrel );
86  clusterAlgo_.setThreshSeedBarrel( threshSeedBarrel );
87 
88  clusterAlgo_.setThreshPtBarrel( threshPtBarrel );
89  clusterAlgo_.setThreshPtSeedBarrel( threshPtSeedBarrel );
90 
91  clusterAlgo_.setThreshCleanBarrel(threshCleanBarrel);
92  clusterAlgo_.setS4S1CleanBarrel(minS4S1CleanBarrel);
93 
94  clusterAlgo_.setThreshDoubleSpikeBarrel( threshDoubleSpikeBarrel );
95  clusterAlgo_.setS6S2DoubleSpikeBarrel( minS6S2DoubleSpikeBarrel );
96 
97  clusterAlgo_.setThreshEndcap( threshEndcap );
98  clusterAlgo_.setThreshSeedEndcap( threshSeedEndcap );
99 
100  clusterAlgo_.setThreshPtEndcap( threshPtEndcap );
101  clusterAlgo_.setThreshPtSeedEndcap( threshPtSeedEndcap );
102 
103  clusterAlgo_.setThreshCleanEndcap(threshCleanEndcap);
104  clusterAlgo_.setS4S1CleanEndcap(minS4S1CleanEndcap);
105 
106  clusterAlgo_.setThreshDoubleSpikeEndcap( threshDoubleSpikeEndcap );
107  clusterAlgo_.setS6S2DoubleSpikeEndcap( minS6S2DoubleSpikeEndcap );
108 
109  clusterAlgo_.setNNeighbours( nNeighbours );
110 
111  // p1 set to the minimum rechit threshold:
112  double posCalcP1 = threshBarrel<threshEndcap ? threshBarrel:threshEndcap;
113  clusterAlgo_.setPosCalcP1( posCalcP1 );
114  clusterAlgo_.setPosCalcNCrystal( posCalcNCrystal );
115  clusterAlgo_.setShowerSigma( showerSigma );
116 
117  clusterAlgo_.setUseCornerCells( useCornerCells );
118  clusterAlgo_.setCleanRBXandHPDs( cleanRBXandHPDs);
119 
120  int dcormode =
121  iConfig.getParameter<int>("depthCor_Mode");
122 
123  double dcora =
124  iConfig.getParameter<double>("depthCor_A");
125  double dcorb =
126  iConfig.getParameter<double>("depthCor_B");
127  double dcorap =
128  iConfig.getParameter<double>("depthCor_A_preshower");
129  double dcorbp =
130  iConfig.getParameter<double>("depthCor_B_preshower");
131 
132  if( dcormode !=0 )
134  dcora, dcorb,
135  dcorap, dcorbp );
136 
137 
138  // access to the collections of rechits from the various detectors:
139 
140 
141  inputTagPFRecHits_ =
142  iConfig.getParameter<InputTag>("PFRecHits");
143  //---ab
144 
145  //inputTagClusterCollectionName_ = iConfig.getParameter<string>("PFClusterCollectionName");
146 
147  // produces<reco::PFClusterCollection>(inputTagClusterCollectionName_);
148  produces<reco::PFClusterCollection>();
149  produces<reco::PFRecHitCollection>("Cleaned");
150 
151  //---ab
152 }
153 
154 
155 
157 
158 
159 
160 
162  const edm::EventSetup& iSetup) {
163 
164 
166 
167  // access the rechits in the event
168  bool found = iEvent.getByLabel( inputTagPFRecHits_, rechitsHandle );
169 
170  if(!found ) {
171 
172  ostringstream err;
173  err<<"cannot find rechits: "<<inputTagPFRecHits_;
174  LogError("PFClusterProducer")<<err.str()<<endl;
175 
176  throw cms::Exception( "MissingProduct", err.str());
177  }
178 
179 
180  // do clustering
181  clusterAlgo_.doClustering( rechitsHandle );
182 
183  if( verbose_ ) {
184  LogInfo("PFClusterProducer")
185  <<" clusters --------------------------------- "<<endl
186  <<clusterAlgo_<<endl;
187  }
188 
189  // get clusters out of the clustering algorithm
190  // and put them in the event. There is no copy.
191  auto_ptr< vector<reco::PFCluster> > outClusters( clusterAlgo_.clusters() );
192  auto_ptr< vector<reco::PFRecHit> > recHitsCleaned ( clusterAlgo_.rechitsCleaned() );
193  iEvent.put( outClusters );
194  iEvent.put( recHitsCleaned, "Cleaned" );
195 
196 }
197 
198 
199 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
static void setDepthCorParameters(int mode, double a, double b, double ap, double bp)
Definition: PFCluster.h:104
virtual void produce(edm::Event &, const edm::EventSetup &)
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
PFClusterProducer(const edm::ParameterSet &)