CMS 3D CMS Logo

SiStripGainFromData.cc
Go to the documentation of this file.
1 // Original Author: Loic QUERTENMONT
2 // Created: Wed Feb 6 08:55:18 CET 2008
3 
4 #include <memory>
5 
14 
20 
27 
29 
32 
45 
48 
51 
54 
57 
58 
59 #include "TFile.h"
60 #include "TObjString.h"
61 #include "TString.h"
62 #include "TH1F.h"
63 #include "TH2F.h"
64 #include "TProfile.h"
65 #include "TF1.h"
66 #include "TROOT.h"
67 
68 #include <unordered_map>
69 #include <memory>
70 
71 
72 
73 using namespace edm;
74 using namespace reco;
75 using namespace std;
76 
77 
78 struct stAPVGain{unsigned int Index; int DetId; int APVId; int SubDet; float Eta; float R; float Phi; float Thickness; double MPV; double Gain; double PreviousGain; char Side;};
79 
80 class SiStripGainFromData : public ConditionDBWriter<SiStripApvGain> {
81  public:
82  explicit SiStripGainFromData(const edm::ParameterSet&);
83  ~SiStripGainFromData() override;
84 
85 
86  private:
87  void algoBeginJob(const edm::EventSetup&) override ;
88  void algoEndJob() override ;
89  void algoBeginRun(const edm::Run &, const edm::EventSetup &) override;
90 // virtual void algoBeginRun(const edm::Event& iEvent, const edm::EventSetup& iSetup);
91  void algoAnalyze(const edm::Event &, const edm::EventSetup &) override;
92 
93  std::unique_ptr<SiStripApvGain> getNewObject() override;
96 
97  double ComputeChargeOverPath(const SiStripCluster* Cluster,TrajectoryStateOnSurface trajState, const edm::EventSetup* iSetup, const Track* track, double trajChi2OverN);
98  bool IsFarFromBorder(TrajectoryStateOnSurface trajState, const uint32_t detid, const edm::EventSetup* iSetup);
99 
100  void getPeakOfLandau(TH1* InputHisto, double* FitResults, double LowRange=0, double HighRange=5400);
101 
104 
106  unsigned int MinNrEntries;
107  double MaxMPVError;
111  double MinTrackEta;
112  double MaxTrackEta;
113  unsigned int MaxNrStrips;
114  unsigned int MinTrackHits;
121 
127 
128  vector<string> VInputFiles;
129 
131 
134 
138  TH1F* HTrackHits;
139 
146 
153 
160 
167 
168 
169  TH1F* Charge_TIB;
170  TH1F* Charge_TID;
171  TH1F* Charge_TIDP;
172  TH1F* Charge_TIDM;
173  TH1F* Charge_TOB;
174  TH1F* Charge_TEC;
175  TH1F* Charge_TEC1;
176  TH1F* Charge_TEC2;
177  TH1F* Charge_TECP;
178  TH1F* Charge_TECM;
179 
180  TH2F* MPV_Vs_Phi;
181  TH2F* MPV_Vs_Eta;
182  TH2F* MPV_Vs_R;
183 
184 // TH2F* PD_Vs_Eta;
185 // TH2F* PD_Vs_R;
186 
187  TH1F* APV_DetId;
188  TH1F* APV_Id;
189  TH1F* APV_Eta;
190  TH1F* APV_R;
191  TH1F* APV_SubDet;
193  TH2F* APV_Charge;
196  TH1F* APV_MPV;
197  TH1F* APV_Gain;
201 
202  TH1F* MPVs;
203  TH1F* MPVs320;
204  TH1F* MPVs500;
205 
206 // TH2F* MPV_vs_10RplusEta;
207 
208 
211 // TH2F* Charge_Vs_PathLength_CS1;
212 // TH2F* Charge_Vs_PathLength_CS2;
213 // TH2F* Charge_Vs_PathLength_CS3;
214 // TH2F* Charge_Vs_PathLength_CS4;
215 // TH2F* Charge_Vs_PathLength_CS5;
216 
217 // TH1F* MPV_Vs_PathLength_CS1;
218 // TH1F* MPV_Vs_PathLength_CS2;
219 // TH1F* MPV_Vs_PathLength_CS3;
220 // TH1F* MPV_Vs_PathLength_CS4;
221 // TH1F* MPV_Vs_PathLength_CS5;
222 
223 // TH1F* FWHM_Vs_PathLength_CS1;
224 // TH1F* FWHM_Vs_PathLength_CS2;
225 // TH1F* FWHM_Vs_PathLength_CS3;
226 // TH1F* FWHM_Vs_PathLength_CS4;
227 // TH1F* FWHM_Vs_PathLength_CS5;
228 
229 
233 
237 
241 
242 
246 
250 
252  TH1F* MPV_Vs_Beta;
254 
259 
261 
264 
265  TH1F* JobInfo;
266 
267  TH1F* HFirstStrip;
268 
269  unsigned int NEvent;
270  unsigned int SRun;
271  unsigned int SEvent;
273  unsigned int ERun;
274  unsigned int EEvent;
276 
277  private :
278  class isEqual{
279  public:
280  template <class T> bool operator () (const T& PseudoDetId1, const T& PseudoDetId2) { return PseudoDetId1==PseudoDetId2; }
281  };
282 
283  std::vector<stAPVGain*> APVsCollOrdered;
284  std::unordered_map<unsigned int, stAPVGain*> APVsColl;
285 };
286 
288 {
289  AlgoMode = iConfig.getParameter<std::string>("AlgoMode");
290 
291  OutputGains = iConfig.getParameter<std::string>("OutputGains");
292  OutputHistos = iConfig.getParameter<std::string>("OutputHistos");
293 
294  TrajToTrackProducer = iConfig.getParameter<std::string>("TrajToTrackProducer");
295  TrajToTrackLabel = iConfig.getParameter<std::string>("TrajToTrackLabel");
296 
297  CheckLocalAngle = iConfig.getUntrackedParameter<bool> ("checkLocalAngle" , false);
298  MinNrEntries = iConfig.getUntrackedParameter<unsigned>("minNrEntries" , 20);
299  MaxMPVError = iConfig.getUntrackedParameter<double> ("maxMPVError" , 500.0);
300  MaxChi2OverNDF = iConfig.getUntrackedParameter<double> ("maxChi2OverNDF" , 5.0);
301  MinTrackMomentum = iConfig.getUntrackedParameter<double> ("minTrackMomentum" , 3.0);
302  MaxTrackMomentum = iConfig.getUntrackedParameter<double> ("maxTrackMomentum" , 99999.0);
303  MinTrackEta = iConfig.getUntrackedParameter<double> ("minTrackEta" , -5.0);
304  MaxTrackEta = iConfig.getUntrackedParameter<double> ("maxTrackEta" , 5.0);
305  MaxNrStrips = iConfig.getUntrackedParameter<unsigned>("maxNrStrips" , 2);
306  MinTrackHits = iConfig.getUntrackedParameter<unsigned>("MinTrackHits" , 8);
307  MaxTrackChiOverNdf = iConfig.getUntrackedParameter<double> ("MaxTrackChiOverNdf" , 3);
308  AllowSaturation = iConfig.getUntrackedParameter<bool> ("AllowSaturation" , false);
309  FirstSetOfConstants = iConfig.getUntrackedParameter<bool> ("FirstSetOfConstants", true);
310  Validation = iConfig.getUntrackedParameter<bool> ("Validation" , false);
311  CheckIfFileExist = iConfig.getUntrackedParameter<bool> ("CheckIfFileExist" , false);
312 
313  CalibrationLevel = iConfig.getUntrackedParameter<int> ("CalibrationLevel" , 0);
314 
315 
316  if( strcmp(AlgoMode.c_str(),"WriteOnDB")==0 )
317  VInputFiles = iConfig.getParameter<vector<string> >("VInputFiles");
318 
321 
322  //if( OutputHistos!="" )
323  // dqmStore_->open(OutputHistos.c_str(), true);
324 }
325 
326 
328 {
329 }
330 
331 
332 void
334 {
335  //Retrieve tracker topology from geometry
336  edm::ESHandle<TrackerTopology> tTopoHandle;
337  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
338  const TrackerTopology* const tTopo = tTopoHandle.product();
339 
340  iSetup_ = &iSetup;
341 
342 // TH1::AddDirectory(kTRUE);
343 
344  tmp = dqmStore_->book1D ("JobInfo" , "JobInfo", 20,0,20); JobInfo = tmp->getTH1F();
345 
346  tmp = dqmStore_->book1D ("APV_DetId" , "APV_DetId" , 72785,0,72784); APV_DetId = tmp->getTH1F();
347  tmp = dqmStore_->book1D ("APV_Id" , "APV_Id" , 72785,0,72784); APV_Id = tmp->getTH1F();
348  tmp = dqmStore_->book1D ("APV_Eta" , "APV_Eta" , 72785,0,72784); APV_Eta = tmp->getTH1F();
349  tmp = dqmStore_->book1D ("APV_R" , "APV_R" , 72785,0,72784); APV_R = tmp->getTH1F();
350  tmp = dqmStore_->book1D ("APV_SubDet" , "APV_SubDet" , 72785,0,72784); APV_SubDet = tmp->getTH1F();
351  tmp = dqmStore_->book2D ("APV_Momentum" , "APV_Momentum" , 72785,0,72784, 50,0,200); APV_Momentum = tmp->getTH2F();
352  tmp = dqmStore_->book2D ("APV_Charge" , "APV_Charge" , 72785,0,72784, 1000,0,2000); APV_Charge = tmp->getTH2F();
353  tmp = dqmStore_->book2D ("APV_PathLength" , "APV_PathLength" , 72785,0,72784, 100,0.2,1.4); APV_PathLength = tmp->getTH2F();
354  tmp = dqmStore_->book1D ("APV_PathLengthM", "APV_PathLengthM", 72785,0,72784); APV_PathLengthM = tmp->getTH1F();
355  tmp = dqmStore_->book1D ("APV_MPV" , "APV_MPV" , 72785,0,72784); APV_MPV = tmp->getTH1F();
356  tmp = dqmStore_->book1D ("APV_Gain" , "APV_Gain" , 72785,0,72784); APV_Gain = tmp->getTH1F();
357  tmp = dqmStore_->book1D ("APV_PrevGain" , "APV_PrevGain" , 72785,0,72784); APV_PrevGain = tmp->getTH1F();
358  tmp = dqmStore_->book1D ("APV_CumulGain" , "APV_CumulGain" , 72785,0,72784); APV_CumulGain = tmp->getTH1F();
359  tmp = dqmStore_->book1D ("APV_Thickness" , "APV_Thicknes" , 72785,0,72784); APV_Thickness = tmp->getTH1F();
360 
361 
362  tmp = dqmStore_->book2D ("Tracks_P_Vs_Eta" , "Tracks_P_Vs_Eta" , 30, 0,3,50,0,200); Tracks_P_Vs_Eta = tmp->getTH2F();
363  tmp = dqmStore_->book2D ("Tracks_Pt_Vs_Eta" , "Tracks_Pt_Vs_Eta", 30, 0,3,50,0,200); Tracks_Pt_Vs_Eta = tmp->getTH2F();
364 
365  tmp = dqmStore_->book2D ("Charge_Vs_PathTIB" , "Charge_Vs_PathTIB" ,100,0.2,1.4, 500,0,2000); Charge_Vs_PathTIB = tmp->getTH2F();
366  tmp = dqmStore_->book2D ("Charge_Vs_PathTID" , "Charge_Vs_PathTID" ,100,0.2,1.4, 500,0,2000); Charge_Vs_PathTID = tmp->getTH2F();
367  tmp = dqmStore_->book2D ("Charge_Vs_PathTOB" , "Charge_Vs_PathTOB" ,100,0.2,1.4, 500,0,2000); Charge_Vs_PathTOB = tmp->getTH2F();
368  tmp = dqmStore_->book2D ("Charge_Vs_PathTEC" , "Charge_Vs_PathTEC" ,100,0.2,1.4, 500,0,2000); Charge_Vs_PathTEC = tmp->getTH2F();
369  tmp = dqmStore_->book2D ("Charge_Vs_PathTEC1", "Charge_Vs_PathTEC1",100,0.2,1.4, 500,0,2000); Charge_Vs_PathTEC1 = tmp->getTH2F();
370  tmp = dqmStore_->book2D ("Charge_Vs_PathTEC2", "Charge_Vs_PathTEC2",100,0.2,1.4, 500,0,2000); Charge_Vs_PathTEC2 = tmp->getTH2F();
371 
372 
373  tmp = dqmStore_->book1D ("Charge_TIB" , "Charge_TIB" ,1000,0,2000); Charge_TIB = tmp->getTH1F();
374  tmp = dqmStore_->book1D ("Charge_TID" , "Charge_TID" ,1000,0,2000); Charge_TID = tmp->getTH1F();
375  tmp = dqmStore_->book1D ("Charge_TID+", "Charge_TID+",1000,0,2000); Charge_TIDP = tmp->getTH1F();
376  tmp = dqmStore_->book1D ("Charge_TID-", "Charge_TID-",1000,0,2000); Charge_TIDM = tmp->getTH1F();
377  tmp = dqmStore_->book1D ("Charge_TOB" , "Charge_TOB" ,1000,0,2000); Charge_TOB = tmp->getTH1F();
378  tmp = dqmStore_->book1D ("Charge_TEC" , "Charge_TEC" ,1000,0,2000); Charge_TEC = tmp->getTH1F();
379  tmp = dqmStore_->book1D ("Charge_TEC1", "Charge_TEC1",1000,0,2000); Charge_TEC1 = tmp->getTH1F();
380  tmp = dqmStore_->book1D ("Charge_TEC2", "Charge_TEC2",1000,0,2000); Charge_TEC2 = tmp->getTH1F();
381  tmp = dqmStore_->book1D ("Charge_TEC+", "Charge_TEC+",1000,0,2000); Charge_TECP = tmp->getTH1F();
382  tmp = dqmStore_->book1D ("Charge_TEC-", "Charge_TEC-",1000,0,2000); Charge_TECM = tmp->getTH1F();
383 
384 
385 /*
386  tmp = dqmStore_->book2D ("Charge_Vs_PathLength_CS1", "Charge_Vs_PathLength_CS1" , 250,0.2,1.4, 500,0,2000); Charge_Vs_PathLength_CS1 = tmp->getTH2F();
387  tmp = dqmStore_->book2D ("Charge_Vs_PathLength_CS2", "Charge_Vs_PathLength_CS2" , 250,0.2,1.4, 500,0,2000); Charge_Vs_PathLength_CS2 = tmp->getTH2F();
388  tmp = dqmStore_->book2D ("Charge_Vs_PathLength_CS3", "Charge_Vs_PathLength_CS3" , 250,0.2,1.4, 500,0,2000); Charge_Vs_PathLength_CS3 = tmp->getTH2F();
389  tmp = dqmStore_->book2D ("Charge_Vs_PathLength_CS4", "Charge_Vs_PathLength_CS4" , 250,0.2,1.4, 500,0,2000); Charge_Vs_PathLength_CS4 = tmp->getTH2F();
390  tmp = dqmStore_->book2D ("Charge_Vs_PathLength_CS5", "Charge_Vs_PathLength_CS5" , 250,0.2,1.4, 500,0,2000); Charge_Vs_PathLength_CS5 = tmp->getTH2F();
391 */
392  tmp = dqmStore_->book2D ("Charge_Vs_PathLength" , "Charge_Vs_PathLength" , 100,0.2,1.4, 500,0,2000); Charge_Vs_PathLength = tmp->getTH2F();
393  tmp = dqmStore_->book2D ("Charge_Vs_PathLength320" , "Charge_Vs_PathLength" , 100,0.2,1.4, 500,0,2000); Charge_Vs_PathLength320 = tmp->getTH2F();
394  tmp = dqmStore_->book2D ("Charge_Vs_PathLength500" , "Charge_Vs_PathLength" , 100,0.2,1.4, 500,0,2000); Charge_Vs_PathLength500 = tmp->getTH2F();
395 
396  tmp = dqmStore_->book2D ("Charge_Vs_TransversAngle" , "Charge_Vs_TransversAngle" , 220,-20,200, 500,0,2000); Charge_Vs_TransversAngle = tmp->getTH2F();
397  tmp = dqmStore_->book2D ("Charge_Vs_Alpha" , "Charge_Vs_Alpha" , 220,-20,200, 500,0,2000); Charge_Vs_Alpha = tmp->getTH2F();
398  tmp = dqmStore_->book2D ("Charge_Vs_Beta" , "Charge_Vs_Beta" , 220,-20,200, 500,0,2000); Charge_Vs_Beta = tmp->getTH2F();
399 
400  tmp = dqmStore_->book2D ("NStrips_Vs_TransversAngle", "NStrips_Vs_TransversAngle", 220,-20,200, 10,0,10); NStrips_Vs_TransversAngle = tmp->getTH2F();
401  tmp = dqmStore_->book2D ("NStrips_Vs_Alpha" , "NStrips_Vs_Alpha" , 220,-20,200, 10,0,10); NStrips_Vs_Alpha = tmp->getTH2F();
402  tmp = dqmStore_->book2D ("NStrips_Vs_Beta" , "NStrips_Vs_Beta" , 220,-20,200, 10,0,10); NStrips_Vs_Beta = tmp->getTH2F();
403  tmp = dqmStore_->book1D ("NHighStripInCluster" , "NHighStripInCluster" , 15,0,14); NHighStripInCluster = tmp->getTH1F();
404  tmp = dqmStore_->book1D ("NSatStripInCluster" , "NSatStripInCluster" , 50,0,50); NSatStripInCluster = tmp->getTH1F();
405 
406  tmp = dqmStore_->book1D ("TrackChi2OverNDF","TrackChi2OverNDF", 500, 0,10); HTrackChi2OverNDF = tmp->getTH1F();
407  tmp = dqmStore_->book1D ("TrackHits","TrackHits", 40, 0,40); HTrackHits = tmp->getTH1F();
408 
409  tmp = dqmStore_->book1D ("FirstStrip","FirstStrip", 800, 0,800); HFirstStrip = tmp->getTH1F();
410 
411  if( strcmp(AlgoMode.c_str(),"MultiJob")!=0 ){
412 
413  tmp = dqmStore_->book2D ("MPV_Vs_EtaTIB" , "MPV_Vs_EtaTIB" , 50, -3.0, 3.0, 300, 0, 600); MPV_Vs_EtaTIB = tmp->getTH2F();
414  tmp = dqmStore_->book2D ("MPV_Vs_EtaTID" , "MPV_Vs_EtaTID" , 50, -3.0, 3.0, 300, 0, 600); MPV_Vs_EtaTID = tmp->getTH2F();
415  tmp = dqmStore_->book2D ("MPV_Vs_EtaTOB" , "MPV_Vs_EtaTOB" , 50, -3.0, 3.0, 300, 0, 600); MPV_Vs_EtaTOB = tmp->getTH2F();
416  tmp = dqmStore_->book2D ("MPV_Vs_EtaTEC" , "MPV_Vs_EtaTEC" , 50, -3.0, 3.0, 300, 0, 600); MPV_Vs_EtaTEC = tmp->getTH2F();
417  tmp = dqmStore_->book2D ("MPV_Vs_EtaTEC1" , "MPV_Vs_EtaTEC1", 50, -3.0, 3.0, 300, 0, 600); MPV_Vs_EtaTEC1 = tmp->getTH2F();
418  tmp = dqmStore_->book2D ("MPV_Vs_EtaTEC2" , "MPV_Vs_EtaTEC2", 50, -3.0, 3.0, 300, 0, 600); MPV_Vs_EtaTEC2 = tmp->getTH2F();
419 
420  tmp = dqmStore_->book2D ("MPV_Vs_PhiTIB" , "MPV_Vs_PhiTIB" , 50, -3.2, 3.2, 300, 0, 600); MPV_Vs_PhiTIB = tmp->getTH2F();
421  tmp = dqmStore_->book2D ("MPV_Vs_PhiTID" , "MPV_Vs_PhiTID" , 50, -3.2, 3.2, 300, 0, 600); MPV_Vs_PhiTID = tmp->getTH2F();
422  tmp = dqmStore_->book2D ("MPV_Vs_PhiTOB" , "MPV_Vs_PhiTOB" , 50, -3.2, 3.2, 300, 0, 600); MPV_Vs_PhiTOB = tmp->getTH2F();
423  tmp = dqmStore_->book2D ("MPV_Vs_PhiTEC" , "MPV_Vs_PhiTEC" , 50, -3.2, 3.2, 300, 0, 600); MPV_Vs_PhiTEC = tmp->getTH2F();
424  tmp = dqmStore_->book2D ("MPV_Vs_PhiTEC1" , "MPV_Vs_PhiTEC1", 50, -3.2, 3.2, 300, 0, 600); MPV_Vs_PhiTEC1 = tmp->getTH2F();
425  tmp = dqmStore_->book2D ("MPV_Vs_PhiTEC2" , "MPV_Vs_PhiTEC2", 50, -3.2, 3.2, 300, 0, 600); MPV_Vs_PhiTEC2 = tmp->getTH2F();
426 
427 
428  tmp = dqmStore_->book1D ("MPV_Vs_PathTIB" , "MPV_Vs_PathTIB" ,100,0.2,1.4); MPV_Vs_PathTIB = tmp->getTH1F();
429  tmp = dqmStore_->book1D ("MPV_Vs_PathTID" , "MPV_Vs_PathTID" ,100,0.2,1.4); MPV_Vs_PathTID = tmp->getTH1F();
430  tmp = dqmStore_->book1D ("MPV_Vs_PathTOB" , "MPV_Vs_PathTOB" ,100,0.2,1.4); MPV_Vs_PathTOB = tmp->getTH1F();
431  tmp = dqmStore_->book1D ("MPV_Vs_PathTEC" , "MPV_Vs_PathTEC" ,100,0.2,1.4); MPV_Vs_PathTEC = tmp->getTH1F();
432  tmp = dqmStore_->book1D ("MPV_Vs_PathTEC1" , "MPV_Vs_PathTEC1" ,100,0.2,1.4); MPV_Vs_PathTEC1 = tmp->getTH1F();
433  tmp = dqmStore_->book1D ("MPV_Vs_PathTEC2" , "MPV_Vs_PathTEC2" ,100,0.2,1.4); MPV_Vs_PathTEC2 = tmp->getTH1F();
434 
435  tmp = dqmStore_->book2D ("MPV_Vs_Phi", "MPV_Vs_Phi", 50, -3.2, 3.2 , 300, 0, 600); MPV_Vs_Phi = tmp->getTH2F();
436  tmp = dqmStore_->book2D ("MPV_Vs_Eta", "MPV_Vs_Eta", 50, -3.0, 3.0 , 300, 0, 600); MPV_Vs_Eta = tmp->getTH2F();
437  tmp = dqmStore_->book2D ("MPV_Vs_R" , "MPV_Vs_R" , 150, 0.0, 150.0, 300, 0, 600); MPV_Vs_R = tmp->getTH2F();
438 /*
439  tmp = dqmStore_->book1D ("MPV_Vs_PathLength_CS1" , "MPV_Vs_PathLength_CS1" , 250, 0.2, 1.4); MPV_Vs_PathLength_CS1 = tmp->getTH1F();
440  tmp = dqmStore_->book1D ("MPV_Vs_PathLength_CS2" , "MPV_Vs_PathLength_CS2" , 250, 0.2, 1.4); MPV_Vs_PathLength_CS2 = tmp->getTH1F();
441  tmp = dqmStore_->book1D ("MPV_Vs_PathLength_CS3" , "MPV_Vs_PathLength_CS3" , 250, 0.2, 1.4); MPV_Vs_PathLength_CS3 = tmp->getTH1F();
442  tmp = dqmStore_->book1D ("MPV_Vs_PathLength_CS4" , "MPV_Vs_PathLength_CS4" , 250, 0.2, 1.4); MPV_Vs_PathLength_CS4 = tmp->getTH1F();
443  tmp = dqmStore_->book1D ("MPV_Vs_PathLength_CS5" , "MPV_Vs_PathLength_CS5" , 250, 0.2, 1.4); MPV_Vs_PathLength_CS5 = tmp->getTH1F();
444 
445  tmp = dqmStore_->book1D ("FWHM_Vs_PathLength_CS1" , "FWHM_Vs_PathLength_CS1", 250, 0.2, 1.4); FWHM_Vs_PathLength_CS1 = tmp->getTH1F();
446  tmp = dqmStore_->book1D ("FWHM_Vs_PathLength_CS2" , "FWHM_Vs_PathLength_CS2", 250, 0.2, 1.4); FWHM_Vs_PathLength_CS2 = tmp->getTH1F();
447  tmp = dqmStore_->book1D ("FWHM_Vs_PathLength_CS3" , "FWHM_Vs_PathLength_CS3", 250, 0.2, 1.4); FWHM_Vs_PathLength_CS3 = tmp->getTH1F();
448  tmp = dqmStore_->book1D ("FWHM_Vs_PathLength_CS4" , "FWHM_Vs_PathLength_CS4", 250, 0.2, 1.4); FWHM_Vs_PathLength_CS4 = tmp->getTH1F();
449  tmp = dqmStore_->book1D ("FWHM_Vs_PathLength_CS5" , "FWHM_Vs_PathLength_CS5", 250, 0.2, 1.4); FWHM_Vs_PathLength_CS5 = tmp->getTH1F();
450 */
451  tmp = dqmStore_->book1D ("MPV_Vs_PathLength" , "MPV_Vs_PathLength" , 100, 0.2, 1.4); MPV_Vs_PathLength = tmp->getTH1F();
452  tmp = dqmStore_->book1D ("MPV_Vs_PathLength320" , "MPV_Vs_PathLength" , 100, 0.2, 1.4); MPV_Vs_PathLength320 = tmp->getTH1F();
453  tmp = dqmStore_->book1D ("MPV_Vs_PathLength500" , "MPV_Vs_PathLength" , 100, 0.2, 1.4); MPV_Vs_PathLength500 = tmp->getTH1F();
454 
455  tmp = dqmStore_->book1D ("FWHM_Vs_PathLength" , "FWHM_Vs_PathLength" , 100, 0.2, 1.4); FWHM_Vs_PathLength = tmp->getTH1F();
456  tmp = dqmStore_->book1D ("FWHM_Vs_PathLength320" , "FWHM_Vs_PathLength" , 100, 0.2, 1.4); FWHM_Vs_PathLength320 = tmp->getTH1F();
457  tmp = dqmStore_->book1D ("FWHM_Vs_PathLength500" , "FWHM_Vs_PathLength" , 100, 0.2, 1.4); FWHM_Vs_PathLength500 = tmp->getTH1F();
458 
459  tmp = dqmStore_->book1D ("MPV_Vs_TransversAngle" , "MPV_Vs_TransversAngle" , 220, -20, 200); MPV_Vs_TransversAngle = tmp->getTH1F();
460  tmp = dqmStore_->book1D ("MPV_Vs_Alpha" , "MPV_Vs_Alpha" , 220, -20, 200); MPV_Vs_Alpha = tmp->getTH1F();
461  tmp = dqmStore_->book1D ("MPV_Vs_Beta" , "MPV_Vs_Beta" , 220, -20, 200); MPV_Vs_Beta = tmp->getTH1F();
462 
463  tmp = dqmStore_->book2D ("Error_Vs_MPV" , "Error_Vs_MPV" ,600,0,600 ,100 ,0 ,50); Error_Vs_MPV = tmp->getTH2F();
464  tmp = dqmStore_->book2D ("Error_Vs_Entries","Error_Vs_Entries",500,0,10000 ,100 ,0 ,50); Error_Vs_Entries = tmp->getTH2F();
465  tmp = dqmStore_->book2D ("Error_Vs_Eta" , "Error_Vs_Eta" ,50 ,-3.0,3.0 ,100 ,0 ,50 ); Error_Vs_Eta = tmp->getTH2F();
466  tmp = dqmStore_->book2D ("Error_Vs_Phi" , "Error_Vs_Phi" ,50 ,-3.2,3.2 ,100 ,0 ,50); Error_Vs_Phi = tmp->getTH2F();
467 
468 
469  tmp = dqmStore_->book2D ("NoMPV_Vs_EtaPhi" , "NoMPV_Vs_EtaPhi" ,50,-3.0,3.0 ,50 ,-3.2,3.2); NoMPV_Vs_EtaPhi = tmp->getTH2F();
470 
471 
472 
473  tmp = dqmStore_->book1D ("NumberOfEntriesByAPV" , "NumberOfEntriesByAPV" , 1000, 0,10000); NumberOfEntriesByAPV = tmp->getTH1F();
474  tmp = dqmStore_->book1D ("Chi2OverNDF","Chi2OverNDF", 500, 0,25); HChi2OverNDF = tmp->getTH1F();
475 
476  tmp = dqmStore_->book1D ("MPVs", "MPVs", 600,0,600); MPVs = tmp->getTH1F();
477  tmp = dqmStore_->book1D ("MPVs320", "MPVs320", 600,0,600); MPVs320 = tmp->getTH1F();
478  tmp = dqmStore_->book1D ("MPVs500", "MPVs500", 600,0,600); MPVs500 = tmp->getTH1F();
479 
480 // MPV_vs_10RplusEta tmp = dqmStore_->book2D ("MPV_vs_10RplusEta","MPV_vs_10RplusEta", 48000,0,2400, 800,100,500);
481  }
482 
483  gROOT->cd();
484 
486  iSetup.get<TrackerDigiGeometryRecord>().get( tkGeom );
487  auto const & Det = tkGeom->dets();
488 
489 
490  edm::ESHandle<SiStripGain> gainHandle;
491 // if(strcmp(AlgoMode.c_str(),"MultiJob")!=0 && !FirstSetOfConstants){
492  iSetup.get<SiStripGainRcd>().get(gainHandle);
493  if(!gainHandle.isValid()){printf("\n#####################\n\nERROR --> gainHandle is not valid\n\n#####################\n\n");exit(0);}
494 // }
495 
496  unsigned int Id=0;
497  for(unsigned int i=0;i<Det.size();i++){
498  DetId Detid = Det[i]->geographicalId();
499  int SubDet = Detid.subdetId();
500 
501  if( SubDet == StripSubdetector::TIB || SubDet == StripSubdetector::TID ||
502  SubDet == StripSubdetector::TOB || SubDet == StripSubdetector::TEC ){
503 
504  auto DetUnit = dynamic_cast<StripGeomDetUnit const*> (Det[i]);
505  if(!DetUnit)continue;
506 
507  const StripTopology& Topo = DetUnit->specificTopology();
508  unsigned int NAPV = Topo.nstrips()/128;
509 
510  double Phi = DetUnit->position().basicVector().phi();
511  double Eta = DetUnit->position().basicVector().eta();
512  double R = DetUnit->position().basicVector().transverse();
513  double Thick = DetUnit->surface().bounds().thickness();
514 
515  for(unsigned int j=0;j<NAPV;j++){
516  stAPVGain* APV = new stAPVGain;
517  APV->Index = Id;
518  APV->DetId = Detid.rawId();
519  APV->APVId = j;
520  APV->SubDet = SubDet;
521  APV->MPV = -1;
522  APV->Gain = -1;
523  APV->PreviousGain = 1;
524  APV->Eta = Eta;
525  APV->Phi = Phi;
526  APV->R = R;
527  APV->Thickness = Thick;
528  APV->Side = 0;
529 
530  if(SubDet==StripSubdetector::TID){
531 
532  APV->Side = tTopo->tecSide(Detid);
533  }else if(SubDet==StripSubdetector::TEC){
534 
535  APV->Side = tTopo->tecSide(Detid);
536  }
537 
538  APVsCollOrdered.push_back(APV);
539  APVsColl[(APV->DetId<<3) | APV->APVId] = APV;
540  Id++;
541 
542  APV_DetId ->Fill(Id,APV->DetId);
543  APV_Id ->Fill(Id,APV->APVId);
544  APV_Eta ->Fill(Id,APV->Eta);
545  APV_R ->Fill(Id,APV->R);
546  APV_SubDet ->Fill(Id,APV->SubDet);
547  APV_Thickness->Fill(Id,APV->Thickness);
548  }
549  }
550  }
551 
552  NEvent = 0;
553  SRun = 0;
554  SEvent = 0;
555  STimestamp = 0;
556  ERun = 0;
557  EEvent = 0;
558  ETimestamp = 0;
559 }
560 
561 void
563 
564 
565  unsigned int I=0;
566 
567  if( strcmp(AlgoMode.c_str(),"WriteOnDB")==0 || strcmp(AlgoMode.c_str(),"Merge")==0){
568  TH1::AddDirectory(kTRUE);
569 
570  TFile* file = nullptr;
571  for(unsigned int f=0;f<VInputFiles.size();f++){
572  printf("Loading New Input File : %s\n", VInputFiles[f].c_str());
573  if(CheckIfFileExist){
574  FILE* doesFileExist = fopen( VInputFiles[f].c_str(), "r" );
575  if(!doesFileExist){
576  printf("File %s doesn't exist\n",VInputFiles[f].c_str());
577  continue;
578  }else{
579  fclose(doesFileExist);
580  }
581  }
582  file = new TFile( VInputFiles[f].c_str() ); if(!file || file->IsZombie() ){printf("### Bug With File %s\n### File will be skipped \n",VInputFiles[f].c_str()); continue;}
583  APV_Charge ->Add( (TH1*) file->FindObjectAny("APV_Charge") , 1);
584  APV_Momentum ->Add( (TH1*) file->FindObjectAny("APV_Momentum") , 1);
585  APV_PathLength ->Add( (TH1*) file->FindObjectAny("APV_PathLength") , 1);
586 
587  Tracks_P_Vs_Eta ->Add( (TH1*) file->FindObjectAny("Tracks_P_Vs_Eta") , 1);
588  Tracks_Pt_Vs_Eta ->Add( (TH1*) file->FindObjectAny("Tracks_Pt_Vs_Eta") , 1);
589 
590  Charge_Vs_PathTIB ->Add( (TH1*) file->FindObjectAny("Charge_Vs_PathTIB") , 1);
591  Charge_Vs_PathTID ->Add( (TH1*) file->FindObjectAny("Charge_Vs_PathTID") , 1);
592  Charge_Vs_PathTOB ->Add( (TH1*) file->FindObjectAny("Charge_Vs_PathTOB") , 1);
593  Charge_Vs_PathTEC ->Add( (TH1*) file->FindObjectAny("Charge_Vs_PathTEC") , 1);
594  Charge_Vs_PathTEC1 ->Add( (TH1*) file->FindObjectAny("Charge_Vs_PathTEC1") , 1);
595  Charge_Vs_PathTEC2 ->Add( (TH1*) file->FindObjectAny("Charge_Vs_PathTEC2") , 1);
596 
597  HTrackChi2OverNDF ->Add( (TH1*) file->FindObjectAny("TrackChi2OverNDF") , 1);
598  HTrackHits ->Add( (TH1*) file->FindObjectAny("TrackHits") , 1);
599 
600  NHighStripInCluster ->Add( (TH1*) file->FindObjectAny("NHighStripInCluster") , 1);
601  NSatStripInCluster ->Add( (TH1*) file->FindObjectAny("NSatStripInCluster") , 1);
602  Charge_Vs_PathLength ->Add( (TH1*) file->FindObjectAny("Charge_Vs_PathLength") , 1);
603  Charge_Vs_PathLength320 ->Add( (TH1*) file->FindObjectAny("Charge_Vs_PathLength320") , 1);
604  Charge_Vs_PathLength500 ->Add( (TH1*) file->FindObjectAny("Charge_Vs_PathLength500") , 1);
605  Charge_Vs_TransversAngle ->Add( (TH1*) file->FindObjectAny("Charge_Vs_TransversAngle") , 1);
606  NStrips_Vs_TransversAngle->Add( (TH1*) file->FindObjectAny("NStrips_Vs_TransversAngle"), 1);
607  Charge_Vs_Alpha ->Add( (TH1*) file->FindObjectAny("Charge_Vs_Alpha") , 1);
608  NStrips_Vs_Alpha ->Add( (TH1*) file->FindObjectAny("NStrips_Vs_Alpha") , 1);
609  HFirstStrip ->Add( (TH1*) file->FindObjectAny("FirstStrip") , 1);
610 
611  TH1F* JobInfo_tmp = (TH1F*) file->FindObjectAny("JobInfo");
612  NEvent += (unsigned int) JobInfo_tmp->GetBinContent(JobInfo_tmp->GetXaxis()->FindBin(1));
613  unsigned int tmp_SRun = (unsigned int) JobInfo_tmp->GetBinContent(JobInfo_tmp->GetXaxis()->FindBin(3));
614  unsigned int tmp_SEvent = (unsigned int) JobInfo_tmp->GetBinContent(JobInfo_tmp->GetXaxis()->FindBin(4));
615  unsigned int tmp_ERun = (unsigned int) JobInfo_tmp->GetBinContent(JobInfo_tmp->GetXaxis()->FindBin(6));
616  unsigned int tmp_EEvent = (unsigned int) JobInfo_tmp->GetBinContent(JobInfo_tmp->GetXaxis()->FindBin(7));
617 
618  if(SRun==0)SRun = tmp_SRun;
619 
620  if (tmp_SRun< SRun){SRun=tmp_SRun; SEvent=tmp_SEvent;}
621  else if(tmp_SRun==SRun && tmp_SEvent<SEvent){SEvent=tmp_SEvent;}
622 
623  if (tmp_ERun> ERun){ERun=tmp_ERun; EEvent=tmp_EEvent;}
624  else if(tmp_ERun==ERun && tmp_EEvent>EEvent){EEvent=tmp_EEvent;}
625 
626  printf("Deleting Current Input File\n");
627  file->Close();
628  delete file;
629  }
630  }
631 
632  JobInfo->Fill(1,NEvent);
633  JobInfo->Fill(3,SRun);
634  JobInfo->Fill(4,SEvent);
635  JobInfo->Fill(6,ERun);
636  JobInfo->Fill(7,EEvent);
637 
638 
639  if( strcmp(AlgoMode.c_str(),"MultiJob")!=0 ){
640  TH1D* Proj = nullptr;
641  double FitResults[5];
642  I=0;
643  for(auto it = APVsColl.begin();it!=APVsColl.end();it++){
644  if( I%3650==0 ) printf("Fitting Histograms \t %6.2f%%\n",(100.0*I)/APVsColl.size());
645  I++;
646  stAPVGain* APV = it->second;
647 
648  int bin = APV_Charge->GetXaxis()->FindBin(APV->Index);
649  Proj = APV_Charge->ProjectionY(" ",bin,bin,"e");
650  Proj = (TH1D*)Proj->Clone();
651  if(Proj==nullptr)continue;
652 
653  // ADD PROJECTTIONS COMMING FROM THE SECOND APV IN THE PAIR
654  if(CalibrationLevel==1){
655  int SecondAPVId = APV->APVId;
656  if(SecondAPVId%2==0){
657  SecondAPVId = SecondAPVId+1;
658  }else{
659  SecondAPVId = SecondAPVId-1;
660  }
661  stAPVGain* APV2 = APVsColl[(APV->DetId<<3) | SecondAPVId];
662 
663  int bin2 = APV_Charge->GetXaxis()->FindBin(APV2->Index);
664  TH1D* Proj2 = APV_Charge->ProjectionY(" ",bin2,bin2,"e");
665  if(Proj2!=nullptr){
666  Proj->Add(Proj2,1);
667  }
668  }else if(CalibrationLevel>1){
669 // printf("%8i %i--> %4.0f + %4.0f\n",APV->DetId, APV->APVId, 0.0, Proj->GetEntries());
670  for(auto it2 = APVsColl.begin();it2!=APVsColl.end();it2++){
671  stAPVGain* APV2 = it2->second;
672 
673  if(APV2->DetId != APV->DetId)continue;
674  if(APV2->APVId == APV->APVId)continue;
675 
676  int bin2 = APV_Charge->GetXaxis()->FindBin(APV2->Index);
677  TH1D* Proj2 = APV_Charge->ProjectionY(" ",bin2,bin2,"e");
678  if(Proj2!=nullptr){
679 // printf("%8i %i--> %4.0f + %4.0f\n",APV2->DetId, APV2->APVId, Proj->GetEntries(), Proj2->GetEntries());
680  Proj->Add(Proj2,1);
681  }
682  }
683 // printf("%8i %i--> %4.0f Full\n",APV->DetId, APV->APVId, Proj->GetEntries());
684  }
685 
686 
687  //std::cout << "Proj->GetEntries(): " << Proj->GetEntries() << ", Proj->GetMean(): " << Proj->GetMean() << std::endl;
688 
689  getPeakOfLandau(Proj,FitResults);
690  APV->MPV = FitResults[0];
691 // printf("MPV = %f - %f\n",FitResults[0], FitResults[1]);
692  if(FitResults[0]!=-0.5 && FitResults[1]<MaxMPVError){
693  APV_MPV->Fill(APV->Index,APV->MPV);
694  MPVs ->Fill(APV->MPV);
695  if(APV->Thickness<0.04) MPVs320->Fill(APV->MPV);
696  if(APV->Thickness>0.04) MPVs500->Fill(APV->MPV);
697 
698  MPV_Vs_R->Fill(APV->R,APV->MPV);
699  MPV_Vs_Eta->Fill(APV->Eta,APV->MPV);
700  if(APV->SubDet==StripSubdetector::TIB) MPV_Vs_EtaTIB ->Fill(APV->Eta,APV->MPV);
701  if(APV->SubDet==StripSubdetector::TID) MPV_Vs_EtaTID ->Fill(APV->Eta,APV->MPV);
702  if(APV->SubDet==StripSubdetector::TOB) MPV_Vs_EtaTOB ->Fill(APV->Eta,APV->MPV);
703  if(APV->SubDet==StripSubdetector::TEC){ MPV_Vs_EtaTEC ->Fill(APV->Eta,APV->MPV);
704  if(APV->Thickness<0.04) MPV_Vs_EtaTEC1->Fill(APV->Eta,APV->MPV);
705  if(APV->Thickness>0.04) MPV_Vs_EtaTEC2->Fill(APV->Eta,APV->MPV);
706  }
707  MPV_Vs_Phi->Fill(APV->Phi,APV->MPV);
708  if(APV->SubDet==StripSubdetector::TIB) MPV_Vs_PhiTIB ->Fill(APV->Phi,APV->MPV);
709  if(APV->SubDet==StripSubdetector::TID) MPV_Vs_PhiTID ->Fill(APV->Phi,APV->MPV);
710  if(APV->SubDet==StripSubdetector::TOB) MPV_Vs_PhiTOB ->Fill(APV->Phi,APV->MPV);
711  if(APV->SubDet==StripSubdetector::TEC){ MPV_Vs_PhiTEC ->Fill(APV->Phi,APV->MPV);
712  if(APV->Thickness<0.04) MPV_Vs_PhiTEC1->Fill(APV->Phi,APV->MPV);
713  if(APV->Thickness>0.04) MPV_Vs_PhiTEC2->Fill(APV->Phi,APV->MPV);
714  }
715 
716  if(APV->SubDet==StripSubdetector::TIB) Charge_TIB ->Add(Proj,1);
717  if(APV->SubDet==StripSubdetector::TID){ Charge_TID ->Add(Proj,1);
718  if(APV->Side==1) Charge_TIDM->Add(Proj,1);
719  if(APV->Side==2) Charge_TIDP->Add(Proj,1);
720  }
721  if(APV->SubDet==StripSubdetector::TOB) Charge_TOB ->Add(Proj,1);
722  if(APV->SubDet==StripSubdetector::TEC){ Charge_TEC ->Add(Proj,1);
723  if(APV->Thickness<0.04) Charge_TEC1->Add(Proj,1);
724  if(APV->Thickness>0.04) Charge_TEC2->Add(Proj,1);
725  if(APV->Side==1) Charge_TECM->Add(Proj,1);
726  if(APV->Side==2) Charge_TECP->Add(Proj,1);
727  }
728  }
729 
730  if(APV->SubDet==StripSubdetector::TIB) Charge_TIB ->Add(Proj,1);
731  if(APV->SubDet==StripSubdetector::TID){ Charge_TID ->Add(Proj,1);
732  if(APV->Side==1) Charge_TIDM->Add(Proj,1);
733  if(APV->Side==2) Charge_TIDP->Add(Proj,1);
734  }
735  if(APV->SubDet==StripSubdetector::TOB) Charge_TOB ->Add(Proj,1);
736  if(APV->SubDet==StripSubdetector::TEC){ Charge_TEC ->Add(Proj,1);
737  if(APV->Thickness<0.04) Charge_TEC1->Add(Proj,1);
738  if(APV->Thickness>0.04) Charge_TEC2->Add(Proj,1);
739  if(APV->Side==1) Charge_TECM->Add(Proj,1);
740  if(APV->Side==2) Charge_TECP->Add(Proj,1);
741  }
742 
743  if(FitResults[0]!=-0.5){
744  HChi2OverNDF->Fill(FitResults[4]);
745  Error_Vs_MPV->Fill(FitResults[0],FitResults[1]);
746  Error_Vs_Entries->Fill(Proj->GetEntries(),FitResults[1]);
747  Error_Vs_Eta->Fill(APV->Eta,FitResults[1]);
748  Error_Vs_Phi->Fill(APV->Phi,FitResults[1]);
749  }
750  NumberOfEntriesByAPV->Fill(Proj->GetEntries());
751  delete Proj;
752 
753 
754  Proj = APV_PathLength->ProjectionY(" ",bin,bin,"e");
755  if(Proj==nullptr)continue;
756 
757  APV_PathLengthM->SetBinContent(APV->Index, Proj->GetMean(1) );
758  APV_PathLengthM->SetBinError (APV->Index, Proj->GetMeanError(1) );
759 // delete Proj;
760  }
761 
762  unsigned int GOOD = 0;
763  unsigned int BAD = 0;
764  double MPVmean = MPVs->GetMean();
765  MPVmean = 300;
766  for(auto it = APVsColl.begin();it!=APVsColl.end();it++){
767 
768  stAPVGain* APV = it->second;
769  if(APV->MPV>0){
770  APV->Gain = APV->MPV / MPVmean; // APV->MPV;
771  GOOD++;
772  }else{
773  NoMPV_Vs_EtaPhi->Fill(APV->Eta, APV->Phi);
774  APV->Gain = 1;
775  BAD++;
776  }
777  if(APV->Gain<=0) APV->Gain = 1;
778  APV_Gain->Fill(APV->Index,APV->Gain);
779 
780  if(!FirstSetOfConstants) APV->Gain *= APV->PreviousGain;
781  APV_CumulGain->Fill(APV->Index,APV->Gain);
782  }
783 
784  for(int j=0;j<Charge_Vs_PathLength->GetXaxis()->GetNbins();j++){
785  Proj = Charge_Vs_PathLength->ProjectionY(" ",j,j,"e");
786  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
787  MPV_Vs_PathLength->SetBinContent (j, FitResults[0]/Charge_Vs_PathLength->GetXaxis()->GetBinCenter(j));
788  MPV_Vs_PathLength->SetBinError (j, FitResults[1]/Charge_Vs_PathLength->GetXaxis()->GetBinCenter(j));
789  FWHM_Vs_PathLength->SetBinContent(j, FitResults[2]/(FitResults[0]/Charge_Vs_PathLength->GetXaxis()->GetBinCenter(j)) );
790  FWHM_Vs_PathLength->SetBinError (j, FitResults[3]/(FitResults[0]/Charge_Vs_PathLength->GetXaxis()->GetBinCenter(j)) );
791  delete Proj;
792  }
793 
794  for(int j=0;j<Charge_Vs_PathLength320->GetXaxis()->GetNbins();j++){
795  Proj = Charge_Vs_PathLength320->ProjectionY(" ",j,j,"e");
796  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
797  MPV_Vs_PathLength320->SetBinContent (j, FitResults[0]/Charge_Vs_PathLength320->GetXaxis()->GetBinCenter(j));
798  MPV_Vs_PathLength320->SetBinError (j, FitResults[1]/Charge_Vs_PathLength320->GetXaxis()->GetBinCenter(j));
799  FWHM_Vs_PathLength320->SetBinContent(j, FitResults[2]/(FitResults[0]/Charge_Vs_PathLength320->GetXaxis()->GetBinCenter(j)) );
800  FWHM_Vs_PathLength320->SetBinError (j, FitResults[3]/(FitResults[0]/Charge_Vs_PathLength320->GetXaxis()->GetBinCenter(j)) );
801  delete Proj;
802  }
803 
804  for(int j=0;j<Charge_Vs_PathLength500->GetXaxis()->GetNbins();j++){
805  Proj = Charge_Vs_PathLength500->ProjectionY(" ",j,j,"e");
806  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
807  MPV_Vs_PathLength500->SetBinContent (j, FitResults[0]/Charge_Vs_PathLength500->GetXaxis()->GetBinCenter(j));
808  MPV_Vs_PathLength500->SetBinError (j, FitResults[1]/Charge_Vs_PathLength500->GetXaxis()->GetBinCenter(j));
809  FWHM_Vs_PathLength500->SetBinContent(j, FitResults[2]/(FitResults[0]/Charge_Vs_PathLength500->GetXaxis()->GetBinCenter(j) ));
810  FWHM_Vs_PathLength500->SetBinError (j, FitResults[3]/(FitResults[0]/Charge_Vs_PathLength500->GetXaxis()->GetBinCenter(j) ));
811  delete Proj;
812  }
813 /*
814  for(int j=0;j<Charge_Vs_PathLength_CS1->GetXaxis()->GetNbins();j++){
815  Proj = Charge_Vs_PathLength_CS1->ProjectionY(" ",j,j,"e");
816  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
817  MPV_Vs_PathLength_CS1->SetBinContent (j, FitResults[0]/Charge_Vs_PathLength_CS1->GetXaxis()->GetBinCenter(j));
818  MPV_Vs_PathLength_CS1->SetBinError (j, FitResults[1]/Charge_Vs_PathLength_CS1->GetXaxis()->GetBinCenter(j));
819  FWHM_Vs_PathLength_CS1->SetBinContent(j, FitResults[2]/(FitResults[0]/Charge_Vs_PathLength_CS1->GetXaxis()->GetBinCenter(j) ));
820  FWHM_Vs_PathLength_CS1->SetBinError (j, FitResults[3]/(FitResults[0]/Charge_Vs_PathLength_CS1->GetXaxis()->GetBinCenter(j) ));
821  delete Proj;
822  }
823 
824  for(int j=0;j<Charge_Vs_PathLength_CS2->GetXaxis()->GetNbins();j++){
825  Proj = Charge_Vs_PathLength_CS2->ProjectionY(" ",j,j,"e");
826  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
827  MPV_Vs_PathLength_CS2->SetBinContent (j, FitResults[0]/Charge_Vs_PathLength_CS2->GetXaxis()->GetBinCenter(j));
828  MPV_Vs_PathLength_CS2->SetBinError (j, FitResults[1]/Charge_Vs_PathLength_CS2->GetXaxis()->GetBinCenter(j));
829  FWHM_Vs_PathLength_CS2->SetBinContent(j, FitResults[2]/(FitResults[0]/Charge_Vs_PathLength_CS2->GetXaxis()->GetBinCenter(j) ));
830  FWHM_Vs_PathLength_CS2->SetBinError (j, FitResults[3]/(FitResults[0]/Charge_Vs_PathLength_CS2->GetXaxis()->GetBinCenter(j) ));
831  delete Proj;
832  }
833 
834  for(int j=0;j<Charge_Vs_PathLength_CS3->GetXaxis()->GetNbins();j++){
835  Proj = Charge_Vs_PathLength_CS3->ProjectionY(" ",j,j,"e");
836  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
837  MPV_Vs_PathLength_CS3->SetBinContent (j, FitResults[0]/Charge_Vs_PathLength_CS3->GetXaxis()->GetBinCenter(j));
838  MPV_Vs_PathLength_CS3->SetBinError (j, FitResults[1]/Charge_Vs_PathLength_CS3->GetXaxis()->GetBinCenter(j));
839  FWHM_Vs_PathLength_CS3->SetBinContent(j, FitResults[2]/(FitResults[0]/Charge_Vs_PathLength_CS3->GetXaxis()->GetBinCenter(j) ));
840  FWHM_Vs_PathLength_CS3->SetBinError (j, FitResults[3]/(FitResults[0]/Charge_Vs_PathLength_CS3->GetXaxis()->GetBinCenter(j) ));
841  delete Proj;
842  }
843 
844  for(int j=0;j<Charge_Vs_PathLength_CS4->GetXaxis()->GetNbins();j++){
845  Proj = Charge_Vs_PathLength_CS4->ProjectionY(" ",j,j,"e");
846  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
847  MPV_Vs_PathLength_CS4->SetBinContent (j, FitResults[0]/Charge_Vs_PathLength_CS4->GetXaxis()->GetBinCenter(j));
848  MPV_Vs_PathLength_CS4->SetBinError (j, FitResults[1]/Charge_Vs_PathLength_CS4->GetXaxis()->GetBinCenter(j));
849  FWHM_Vs_PathLength_CS4->SetBinContent(j, FitResults[2]/(FitResults[0]/Charge_Vs_PathLength_CS4->GetXaxis()->GetBinCenter(j) ));
850  FWHM_Vs_PathLength_CS4->SetBinError (j, FitResults[3]/(FitResults[0]/Charge_Vs_PathLength_CS4->GetXaxis()->GetBinCenter(j) ));
851  delete Proj;
852  }
853 
854  for(int j=0;j<Charge_Vs_PathLength_CS5->GetXaxis()->GetNbins();j++){
855  Proj = Charge_Vs_PathLength_CS5->ProjectionY(" ",j,j,"e");
856  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
857  MPV_Vs_PathLength_CS5->SetBinContent (j, FitResults[0]/Charge_Vs_PathLength_CS5->GetXaxis()->GetBinCenter(j));
858  MPV_Vs_PathLength_CS5->SetBinError (j, FitResults[1]/Charge_Vs_PathLength_CS5->GetXaxis()->GetBinCenter(j));
859  FWHM_Vs_PathLength_CS5->SetBinContent(j, FitResults[2]/(FitResults[0]/Charge_Vs_PathLength_CS5->GetXaxis()->GetBinCenter(j) ));
860  FWHM_Vs_PathLength_CS5->SetBinError (j, FitResults[3]/(FitResults[0]/Charge_Vs_PathLength_CS5->GetXaxis()->GetBinCenter(j) ));
861  delete Proj;
862  }
863 */
864 
865 
866  for(int j=0;j<Charge_Vs_PathTIB->GetXaxis()->GetNbins();j++){
867  Proj = Charge_Vs_PathTIB->ProjectionY(" ",j,j,"e");
868  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
869  MPV_Vs_PathTIB->SetBinContent(j, FitResults[0]/Charge_Vs_PathTIB->GetXaxis()->GetBinCenter(j));
870  MPV_Vs_PathTIB->SetBinError (j, FitResults[1]/Charge_Vs_PathTIB->GetXaxis()->GetBinCenter(j));
871  delete Proj;
872  }
873 
874  for(int j=0;j<Charge_Vs_PathTID->GetXaxis()->GetNbins();j++){
875  Proj = Charge_Vs_PathTID->ProjectionY(" ",j,j,"e");
876  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
877  MPV_Vs_PathTID->SetBinContent(j, FitResults[0]/Charge_Vs_PathTID->GetXaxis()->GetBinCenter(j));
878  MPV_Vs_PathTID->SetBinError (j, FitResults[1]/Charge_Vs_PathTID->GetXaxis()->GetBinCenter(j));
879  delete Proj;
880  }
881 
882  for(int j=0;j<Charge_Vs_PathTOB->GetXaxis()->GetNbins();j++){
883  Proj = Charge_Vs_PathTOB->ProjectionY(" ",j,j,"e");
884  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
885  MPV_Vs_PathTOB->SetBinContent(j, FitResults[0]/Charge_Vs_PathTOB->GetXaxis()->GetBinCenter(j));
886  MPV_Vs_PathTOB->SetBinError (j, FitResults[1]/Charge_Vs_PathTOB->GetXaxis()->GetBinCenter(j));
887  delete Proj;
888  }
889 
890  for(int j=0;j<Charge_Vs_PathTEC->GetXaxis()->GetNbins();j++){
891  Proj = Charge_Vs_PathTEC->ProjectionY(" ",j,j,"e");
892  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
893  MPV_Vs_PathTEC->SetBinContent(j, FitResults[0]/Charge_Vs_PathTEC->GetXaxis()->GetBinCenter(j));
894  MPV_Vs_PathTEC->SetBinError (j, FitResults[1]/Charge_Vs_PathTEC->GetXaxis()->GetBinCenter(j));
895  delete Proj;
896  }
897 
898  for(int j=0;j<Charge_Vs_PathTEC1->GetXaxis()->GetNbins();j++){
899  Proj = Charge_Vs_PathTEC1->ProjectionY(" ",j,j,"e");
900  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
901  MPV_Vs_PathTEC1->SetBinContent(j, FitResults[0]/Charge_Vs_PathTEC1->GetXaxis()->GetBinCenter(j));
902  MPV_Vs_PathTEC1->SetBinError (j, FitResults[1]/Charge_Vs_PathTEC1->GetXaxis()->GetBinCenter(j));
903  delete Proj;
904  }
905 
906  for(int j=0;j<Charge_Vs_PathTEC2->GetXaxis()->GetNbins();j++){
907  Proj = Charge_Vs_PathTEC2->ProjectionY(" ",j,j,"e");
908  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
909  MPV_Vs_PathTEC2->SetBinContent(j, FitResults[0]/Charge_Vs_PathTEC2->GetXaxis()->GetBinCenter(j));
910  MPV_Vs_PathTEC2->SetBinError (j, FitResults[1]/Charge_Vs_PathTEC2->GetXaxis()->GetBinCenter(j));
911  delete Proj;
912  }
913 
914 
915  for(int j=1;j<Charge_Vs_TransversAngle->GetXaxis()->GetNbins();j++){
916  Proj = Charge_Vs_TransversAngle->ProjectionY(" ",j,j,"e");
917  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
918  MPV_Vs_TransversAngle->SetBinContent(j, FitResults[0]);
919  MPV_Vs_TransversAngle->SetBinError (j, FitResults[1]);
920  delete Proj;
921  }
922 
923  for(int j=1;j<Charge_Vs_Alpha->GetXaxis()->GetNbins();j++){
924  Proj = Charge_Vs_Alpha->ProjectionY(" ",j,j,"e");
925  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
926  MPV_Vs_Alpha->SetBinContent(j, FitResults[0]);
927  MPV_Vs_Alpha->SetBinError (j, FitResults[1]);
928  delete Proj;
929  }
930 
931  for(int j=1;j<Charge_Vs_Beta->GetXaxis()->GetNbins();j++){
932  Proj = Charge_Vs_Beta->ProjectionY(" ",j,j,"e");
933  getPeakOfLandau(Proj,FitResults); if(FitResults[0] ==-0.5)continue;
934  MPV_Vs_Beta->SetBinContent(j, FitResults[0]);
935  MPV_Vs_Beta->SetBinError (j, FitResults[1]);
936  delete Proj;
937  }
938 
939  FILE* Gains = fopen(OutputGains.c_str(),"w");
940  fprintf(Gains,"NEvents = %i\n",NEvent);
941  fprintf(Gains,"Number of APVs = %lu\n",static_cast<unsigned long>(APVsColl.size()));
942  fprintf(Gains,"GoodFits = %i BadFits = %i ratio = %f\n",GOOD,BAD,(100.0*GOOD)/(GOOD+BAD));
943  for(std::vector<stAPVGain*>::iterator it = APVsCollOrdered.begin();it!=APVsCollOrdered.end();it++){
944  stAPVGain* APV = *it;
945  fprintf(Gains,"%i | %i | PreviousGain = %7.5f NewGain = %7.5f\n", APV->DetId,APV->APVId,APV->PreviousGain,APV->Gain);
946  }
947 
948  std::vector<int> DetIdOfBuggedAPV;
949  fprintf(Gains,"----------------------------------------------------------------------\n");
950  for(std::vector<stAPVGain*>::iterator it = APVsCollOrdered.begin();it!=APVsCollOrdered.end();it++){
951  stAPVGain* APV = *it;
952  if(APV->MPV>0 && APV->MPV<200){
953  bool tmpBug = false;
954  for(unsigned int b=0;b<DetIdOfBuggedAPV.size()&&!tmpBug;b++){if(DetIdOfBuggedAPV[b]==APV->DetId)tmpBug=true;}
955  if(!tmpBug){fprintf(Gains,"%i,\n",APV->DetId);DetIdOfBuggedAPV.push_back(APV->DetId);}
956  }
957  }
958 
959 
960  fclose(Gains);
961 
962 // delete [] FitResults;
963 // delete Proj;
964  }
965 
966  dqmStore_->cd();
967  dqmStore_->save(OutputHistos);
968 
969 }
970 
971 
973 
974  edm::ESHandle<SiStripGain> gainHandle;
975  if((strcmp(AlgoMode.c_str(),"MultiJob")!=0 && !FirstSetOfConstants) || Validation){
976  iSetup.get<SiStripGainRcd>().get(gainHandle);
977  if(!gainHandle.isValid()){printf("\n#####################\n\nERROR --> gainHandle is not valid\n\n#####################\n\n");exit(0);}
978 
979 
980  for(std::vector<stAPVGain*>::iterator it = APVsCollOrdered.begin();it!=APVsCollOrdered.end();it++){
981  stAPVGain* APV = *it;
982 
983  if(gainHandle.isValid()){
984  SiStripApvGain::Range detGainRange = gainHandle->getRange(APV->DetId);
985  APV->PreviousGain = *(detGainRange.first + APV->APVId);
986 // APV_PrevGain->Fill(APV->Index,APV->PreviousGain);
987  APV_PrevGain->SetBinContent(APV_PrevGain->GetXaxis()->FindBin(APV->Index),APV->PreviousGain);
988  if(APV->PreviousGain<0)APV->PreviousGain = 1;
989  }else{
990  printf("GAIN HANDLE IS NOT VALID\n");
991  }
992  }
993  }
994 
995 
996 }
997 
998 
999 
1000 void
1002 {
1003 
1004  if( strcmp(AlgoMode.c_str(),"WriteOnDB")==0 ) return;
1005 
1006  if(NEvent==0){
1007  SRun = iEvent.id().run();
1008  SEvent = iEvent.id().event();
1009  STimestamp = iEvent.time().value();
1010  }
1011  ERun = iEvent.id().run();
1012  EEvent = iEvent.id().event();
1013  ETimestamp = iEvent.time().value();
1014  NEvent++;
1015 
1016  iEvent_ = &iEvent;
1017 
1018  Handle<TrajTrackAssociationCollection> trajTrackAssociationHandle;
1019  iEvent.getByLabel(TrajToTrackProducer, TrajToTrackLabel, trajTrackAssociationHandle);
1020  const TrajTrackAssociationCollection TrajToTrackMap = *trajTrackAssociationHandle.product();
1021 
1022 
1023  for(TrajTrackAssociationCollection::const_iterator it = TrajToTrackMap.begin(); it!=TrajToTrackMap.end(); ++it) {
1024  const Track track = *it->val;
1025  const Trajectory traj = *it->key;
1026 
1027 
1028  if(track.p()<MinTrackMomentum || track.p()>MaxTrackMomentum || track.eta()<MinTrackEta || track.eta()>MaxTrackEta) continue;
1029 
1030  Tracks_Pt_Vs_Eta->Fill(fabs(track.eta()),track.pt());
1031  Tracks_P_Vs_Eta ->Fill(fabs(track.eta()),track.p());
1032 
1033  //BEGIN TO COMPUTE NDOF FOR TRACKS NO IMPLEMENTED BEFORE 200pre3
1034  int ndof =0;
1035  const Trajectory::RecHitContainer transRecHits = traj.recHits();
1036 
1037  for(Trajectory::RecHitContainer::const_iterator rechit = transRecHits.begin(); rechit != transRecHits.end(); ++rechit)
1038  if ((*rechit)->isValid()) ndof += (*rechit)->dimension();
1039  ndof -= 5;
1040  //END TO COMPUTE NDOF FOR TRACKS NO IMPLEMENTED BEFORE 200pre3
1041 
1042  HTrackChi2OverNDF->Fill(traj.chiSquared()/ndof);
1043  if(traj.chiSquared()/ndof>MaxTrackChiOverNdf)continue;
1044 
1045  vector<TrajectoryMeasurement> measurements = traj.measurements();
1046  HTrackHits->Fill(traj.foundHits());
1047  if(traj.foundHits()<(int)MinTrackHits)continue;
1048 /*
1049  //BEGIN TO COMPUTE #MATCHEDRECHIT IN THE TRACK
1050  int NMatchedHit = 0;
1051  for(vector<TrajectoryMeasurement>::const_iterator measurement_it = measurements.begin(); measurement_it!=measurements.end(); measurement_it++){
1052  TrajectoryStateOnSurface trajState = measurement_it->updatedState();
1053  if( !trajState.isValid() ) continue;
1054  const TrackingRecHit* hit = (*measurement_it->recHit()).hit();
1055  const SiStripMatchedRecHit2D* sistripmatchedhit = dynamic_cast<const SiStripMatchedRecHit2D*>(hit);
1056  if(sistripmatchedhit)NMatchedHit++;
1057 // NMatchedHit++;
1058 
1059  }
1060  //END TO COMPUTE #MATCHEDRECHIT IN THE TRACK
1061 
1062  if(NMatchedHit<2){
1063  printf("NOT ENOUGH MATCHED RECHITS : %i\n",NMatchedHit);
1064  continue;
1065  }
1066 */
1067 
1068  for(vector<TrajectoryMeasurement>::const_iterator measurement_it = measurements.begin(); measurement_it!=measurements.end(); measurement_it++){
1069 
1070  TrajectoryStateOnSurface trajState = measurement_it->updatedState();
1071  if( !trajState.isValid() ) continue;
1072 
1073  const TrackingRecHit* hit = (*measurement_it->recHit()).hit();
1074  const SiStripRecHit1D* sistripsimple1dhit = dynamic_cast<const SiStripRecHit1D*>(hit);
1075  const SiStripRecHit2D* sistripsimplehit = dynamic_cast<const SiStripRecHit2D*>(hit);
1076  const SiStripMatchedRecHit2D* sistripmatchedhit = dynamic_cast<const SiStripMatchedRecHit2D*>(hit);
1077 
1078  if(sistripsimplehit){
1079  ComputeChargeOverPath((sistripsimplehit->cluster()).get(), trajState, &iSetup, &track, traj.chiSquared()/ndof);
1080  }else if(sistripmatchedhit){
1081  ComputeChargeOverPath(&sistripmatchedhit->monoCluster(),trajState, &iSetup, &track, traj.chiSquared()/ndof);
1082  ComputeChargeOverPath(&sistripmatchedhit->stereoCluster(),trajState, &iSetup, &track, traj.chiSquared()/ndof);
1083  }else if(sistripsimple1dhit){
1084  ComputeChargeOverPath((sistripsimple1dhit->cluster()).get(), trajState, &iSetup, &track, traj.chiSquared()/ndof);
1085  }else{
1086  }
1087 
1088  }
1089 
1090  }
1091 }
1092 
1093 double
1094 SiStripGainFromData::ComputeChargeOverPath(const SiStripCluster* Cluster ,TrajectoryStateOnSurface trajState, const edm::EventSetup* iSetup, const Track* track, double trajChi2OverN)
1095 {
1096  LocalVector trackDirection = trajState.localDirection();
1097  double cosine = trackDirection.z()/trackDirection.mag();
1098  auto const & Ampls = Cluster->amplitudes();
1099  uint32_t DetId = 0; // is 0 since long time Cluster->geographicalId();
1100  int FirstStrip = Cluster->firstStrip();
1101  int APVId = FirstStrip/128;
1102  stAPVGain* APV = APVsColl[(DetId<<3) | APVId];
1103  int Saturation = 0;
1104  bool Overlaping = false;
1105  int Charge = 0;
1106  unsigned int NHighStrip = 0;
1107 
1108  if(!IsFarFromBorder(trajState,DetId, iSetup))return -1;
1109 
1110 
1111  if(FirstStrip==0 )Overlaping=true;
1112  if(FirstStrip==128 )Overlaping=true;
1113  if(FirstStrip==256 )Overlaping=true;
1114  if(FirstStrip==384 )Overlaping=true;
1115  if(FirstStrip==512 )Overlaping=true;
1116  if(FirstStrip==640 )Overlaping=true;
1117 
1118  if(FirstStrip<=127 && FirstStrip+Ampls.size()>127)Overlaping=true;
1119  if(FirstStrip<=255 && FirstStrip+Ampls.size()>255)Overlaping=true;
1120  if(FirstStrip<=383 && FirstStrip+Ampls.size()>383)Overlaping=true;
1121  if(FirstStrip<=511 && FirstStrip+Ampls.size()>511)Overlaping=true;
1122  if(FirstStrip<=639 && FirstStrip+Ampls.size()>639)Overlaping=true;
1123 
1124  if(FirstStrip+Ampls.size()==127 )Overlaping=true;
1125  if(FirstStrip+Ampls.size()==255 )Overlaping=true;
1126  if(FirstStrip+Ampls.size()==383 )Overlaping=true;
1127  if(FirstStrip+Ampls.size()==511 )Overlaping=true;
1128  if(FirstStrip+Ampls.size()==639 )Overlaping=true;
1129  if(FirstStrip+Ampls.size()==767 )Overlaping=true;
1130  if(Overlaping)return -1;
1131 
1132 
1133 /*
1134  if(FirstStrip==0 )Overlaping=true;
1135  if(FirstStrip<=255 && FirstStrip+Ampls.size()>255)Overlaping=true;
1136  if(FirstStrip<=511 && FirstStrip+Ampls.size()>511)Overlaping=true;
1137  if(FirstStrip+Ampls.size()==511 )Overlaping=true;
1138  if(FirstStrip+Ampls.size()==767 )Overlaping=true;
1139  if(Overlaping)return -1;
1140 */
1141 
1142  for(unsigned int a=0;a<Ampls.size();a++){Charge+=Ampls[a];if(Ampls[a]>=254)Saturation++;if(Ampls[a]>=20)NHighStrip++;}
1143  double path = (10.0*APV->Thickness)/fabs(cosine);
1144  double ClusterChargeOverPath = (double)Charge / path ;
1145 
1146  NSatStripInCluster->Fill(Saturation);
1147 
1148  if(Ampls.size()>MaxNrStrips) return -1;
1149  if(Saturation>0 && !AllowSaturation)return -1;
1150  Charge_Vs_PathLength ->Fill(path,Charge);
1151  if(APV->Thickness<0.04) Charge_Vs_PathLength320->Fill(path,Charge);
1152  if(APV->Thickness>0.04) Charge_Vs_PathLength500->Fill(path,Charge);
1153  if(APV->SubDet==StripSubdetector::TIB) Charge_Vs_PathTIB ->Fill(path,Charge);
1154  if(APV->SubDet==StripSubdetector::TID) Charge_Vs_PathTID ->Fill(path,Charge);
1155  if(APV->SubDet==StripSubdetector::TOB) Charge_Vs_PathTOB ->Fill(path,Charge);
1156  if(APV->SubDet==StripSubdetector::TEC){ Charge_Vs_PathTEC ->Fill(path,Charge);
1157  if(APV->Thickness<0.04) Charge_Vs_PathTEC1 ->Fill(path,Charge);
1158  if(APV->Thickness>0.04) Charge_Vs_PathTEC2 ->Fill(path,Charge); }
1159 
1160  double trans = atan2(trackDirection.y(),trackDirection.x())*(180/3.14159265);
1161  double alpha = acos(trackDirection.x() / sqrt( pow(trackDirection.x(),2) + pow(trackDirection.z(),2) ) ) * (180/3.14159265);
1162  double beta = acos(trackDirection.y() / sqrt( pow(trackDirection.x(),2) + pow(trackDirection.z(),2) ) ) * (180/3.14159265);
1163 
1164  if(path>0.4 && path<0.45){
1165  Charge_Vs_TransversAngle->Fill(trans,Charge/path);
1166  Charge_Vs_Alpha->Fill(alpha,Charge/path);
1167  Charge_Vs_Beta->Fill(beta,Charge/path);
1168  }
1169 
1170  NStrips_Vs_TransversAngle->Fill(trans,Ampls.size());
1171  NStrips_Vs_Alpha ->Fill(alpha,Ampls.size());
1172  NStrips_Vs_Beta ->Fill(beta ,Ampls.size());
1173 
1174  NHighStripInCluster->Fill(NHighStrip);
1175 // if(NHighStrip==1) Charge_Vs_PathLength_CS1->Fill(path, Charge );
1176 // if(NHighStrip==2) Charge_Vs_PathLength_CS2->Fill(path, Charge );
1177 // if(NHighStrip==3) Charge_Vs_PathLength_CS3->Fill(path, Charge );
1178 // if(NHighStrip==4) Charge_Vs_PathLength_CS4->Fill(path, Charge );
1179 // if(NHighStrip==5) Charge_Vs_PathLength_CS5->Fill(path, Charge );
1180 
1181  HFirstStrip ->Fill(FirstStrip);
1182 
1183 
1184  if(Validation){ClusterChargeOverPath=ClusterChargeOverPath/APV->PreviousGain;}
1185 
1186  APV_Charge ->Fill(APV->Index,ClusterChargeOverPath);
1187  APV_Momentum ->Fill(APV->Index,trajState.globalMomentum().mag());
1188  APV_PathLength->Fill(APV->Index,path);
1189 
1190  return ClusterChargeOverPath;
1191 }
1192 
1193 bool SiStripGainFromData::IsFarFromBorder(TrajectoryStateOnSurface trajState, const uint32_t detid, const edm::EventSetup* iSetup)
1194 {
1195  edm::ESHandle<TrackerGeometry> tkGeom; iSetup->get<TrackerDigiGeometryRecord>().get( tkGeom );
1196 
1197  LocalPoint HitLocalPos = trajState.localPosition();
1198  LocalError HitLocalError = trajState.localError().positionError() ;
1199 
1200  const GeomDetUnit* it = tkGeom->idToDetUnit(DetId(detid));
1201  if (dynamic_cast<const StripGeomDetUnit*>(it)==nullptr && dynamic_cast<const PixelGeomDetUnit*>(it)==nullptr) {
1202  std::cout << "this detID doesn't seem to belong to the Tracker" << std::endl;
1203  return false;
1204  }
1205 
1206  const BoundPlane plane = it->surface();
1207  const TrapezoidalPlaneBounds* trapezoidalBounds( dynamic_cast<const TrapezoidalPlaneBounds*>(&(plane.bounds())));
1208  const RectangularPlaneBounds* rectangularBounds( dynamic_cast<const RectangularPlaneBounds*>(&(plane.bounds())));
1209 
1210  double DistFromBorder = 1.0;
1211  //double HalfWidth = it->surface().bounds().width() /2.0;
1212  double HalfLength = it->surface().bounds().length() /2.0;
1213 
1214  if(trapezoidalBounds)
1215  {
1216  std::array<const float, 4> const & parameters = (*trapezoidalBounds).parameters();
1217  HalfLength = parameters[3];
1218  //double t = (HalfLength + HitLocalPos.y()) / (2*HalfLength) ;
1219  //HalfWidth = parameters[0] + (parameters[1]-parameters[0]) * t;
1220  }else if(rectangularBounds){
1221  //HalfWidth = it->surface().bounds().width() /2.0;
1222  HalfLength = it->surface().bounds().length() /2.0;
1223  }else{return false;}
1224 
1225 // if (fabs(HitLocalPos.x())+HitLocalError.xx() >= (HalfWidth - DistFromBorder) ) return false;//Don't think is really necessary
1226  if (fabs(HitLocalPos.y())+HitLocalError.yy() >= (HalfLength - DistFromBorder) ) return false;
1227 
1228  return true;
1229 }
1230 
1231 void SiStripGainFromData::getPeakOfLandau(TH1* InputHisto, double* FitResults, double LowRange, double HighRange)
1232 {
1233  double adcs = -0.5;
1234  double adcs_err = 0.;
1235  double width = -0.5;
1236  double width_err = 0;
1237  double chi2overndf = -0.5;
1238 
1239  double nr_of_entries = InputHisto->GetEntries();
1240 
1241  if( (unsigned int)nr_of_entries < MinNrEntries){
1242  FitResults[0] = adcs;
1243  FitResults[1] = adcs_err;
1244  FitResults[2] = width;
1245  FitResults[3] = width_err;
1246  FitResults[4] = chi2overndf;
1247  return;
1248  }
1249 
1250  // perform fit with standard landau
1251  TF1* MyLandau = new TF1("MyLandau","landau",LowRange, HighRange);
1252  MyLandau->SetParameter("MPV",300);
1253 
1254  InputHisto->Fit(MyLandau,"QR WW");
1255 
1256  // MPV is parameter 1 (0=constant, 1=MPV, 2=Sigma)
1257  adcs = MyLandau->GetParameter("MPV");
1258  adcs_err = MyLandau->GetParError(1);
1259  width = MyLandau->GetParameter(2);
1260  width_err = MyLandau->GetParError(2);
1261  chi2overndf = MyLandau->GetChisquare() / MyLandau->GetNDF();
1262 
1263  // if still wrong, give up
1264  if(adcs<2. || chi2overndf>MaxChi2OverNDF){
1265  adcs = -0.5; adcs_err = 0.;
1266  width = -0.5; width_err = 0;
1267  chi2overndf = -0.5;
1268  }
1269 
1270  FitResults[0] = adcs;
1271  FitResults[1] = adcs_err;
1272  FitResults[2] = width;
1273  FitResults[3] = width_err;
1274  FitResults[4] = chi2overndf;
1275 
1276  delete MyLandau;
1277 }
1278 
1279 
1280 
1281 
1282 std::unique_ptr<SiStripApvGain> SiStripGainFromData::getNewObject()
1283 {
1284  cout << "START getNewObject\n";
1285 
1286 // if( !(strcmp(AlgoMode.c_str(),"WriteOnDB")==0 || strcmp(AlgoMode.c_str(),"SingleJob")==0) )return NULL;
1287  if( !(strcmp(AlgoMode.c_str(),"WriteOnDB")==0 || strcmp(AlgoMode.c_str(),"SingleJob")==0) )
1288  return std::make_unique<SiStripApvGain>();
1289 
1290 
1291  auto obj = std::make_unique<SiStripApvGain>();
1292  std::vector<float>* theSiStripVector = nullptr;
1293  int PreviousDetId = -1;
1294  for(unsigned int a=0;a<APVsCollOrdered.size();a++)
1295  {
1297  if(APV==nullptr){ printf("Bug\n"); continue; }
1298  if(APV->DetId != PreviousDetId){
1299  if(theSiStripVector!=nullptr){
1300  SiStripApvGain::Range range(theSiStripVector->begin(),theSiStripVector->end());
1301  if ( !obj->put(PreviousDetId,range) ) printf("Bug to put detId = %i\n",PreviousDetId);
1302  }
1303 
1304  theSiStripVector = new std::vector<float>;
1305  PreviousDetId = APV->DetId;
1306  }
1307  printf("%i | %i | PreviousGain = %7.5f NewGain = %7.5f\n", APV->DetId,APV->APVId,APV->PreviousGain,APV->Gain);
1308  theSiStripVector->push_back(APV->Gain);
1309 // theSiStripVector->push_back(APV->Gain);
1310  }
1311 
1312  if(theSiStripVector!=nullptr){
1313  SiStripApvGain::Range range(theSiStripVector->begin(),theSiStripVector->end());
1314  if ( !obj->put(PreviousDetId,range) ) printf("Bug to put detId = %i\n",PreviousDetId);
1315  }
1316 
1317  cout << "END getNewObject\n";
1318  return obj;
1319 }
1320 
1321 
1322 
RunNumber_t run() const
Definition: EventID.h:39
ClusterRef cluster() const
double p() const
momentum vector magnitude
Definition: TrackBase.h:648
bool IsFarFromBorder(const TrajectoryStateOnSurface &trajState, const GeomDetUnit *it)
Definition: DeDxTools.cc:311
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
int foundHits() const
Definition: Trajectory.h:225
T getUntrackedParameter(std::string const &, T const &) const
float alpha
Definition: AMPTWrapper.h:95
virtual float length() const =0
void algoEndJob() override
SiStripGainFromData(const edm::ParameterSet &)
const_iterator end() const
last iterator over the map (read only)
TH1F * getTH1F() const
std::vector< stAPVGain * > APVsCollOrdered
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
SiStripCluster const & monoCluster() const
LocalVector localDirection() const
unsigned int APVId
Definition: APVGainStruct.h:11
float Thickness
Definition: APVGainStruct.h:19
bool IsFarFromBorder(TrajectoryStateOnSurface trajState, const uint32_t detid, const edm::EventSetup *iSetup)
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
T y() const
Definition: PV3DBase.h:63
void algoBeginRun(const edm::Run &, const edm::EventSetup &) override
ConstRecHitContainer recHits() const
Definition: Trajectory.h:204
const Bounds & bounds() const
Definition: Surface.h:120
uint16_t firstStrip() const
unsigned int DetId
Definition: APVGainStruct.h:10
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
const edm::EventSetup * iSetup_
LocalError positionError() const
const DetContainer & dets() const override
Returm a vector of all GeomDet (including all GeomDetUnits)
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
DataContainer const & measurements() const
Definition: Trajectory.h:196
vector< string > VInputFiles
int iEvent
Definition: GenABIO.cc:230
double Gain
Definition: APVGainStruct.h:26
T mag() const
Definition: PV3DBase.h:67
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:684
std::unordered_map< unsigned int, stAPVGain * > APVsColl
double PreviousGain
Definition: APVGainStruct.h:28
float yy() const
Definition: LocalError.h:26
unsigned int Index
Definition: APVGainStruct.h:8
T sqrt(T t)
Definition: SSEVec.h:18
double pt() const
track transverse momentum
Definition: TrackBase.h:654
T z() const
Definition: PV3DBase.h:64
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
std::pair< ContainerIterator, ContainerIterator > Range
const std::complex< double > I
Definition: I.h:8
double f[11][100]
ClusterRef cluster() const
double const GOOD
Definition: Constants.h:15
const LocalTrajectoryError & localError() const
void algoBeginJob(const edm::EventSetup &) override
unsigned long long TimeValue_t
Definition: Timestamp.h:28
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:535
TH2F * getTH2F() const
bin
set the eta bin as selection string.
Definition: DetId.h:18
double const BAD
Definition: Constants.h:17
ConstRecHitContainer RecHitContainer
Definition: Trajectory.h:46
T const * product() const
Definition: Handle.h:81
double ComputeChargeOverPath(const SiStripCluster *Cluster, TrajectoryStateOnSurface trajState, const edm::EventSetup *iSetup, const Track *track, double trajChi2OverN)
void algoAnalyze(const edm::Event &, const edm::EventSetup &) override
virtual int nstrips() const =0
std::unique_ptr< SiStripApvGain > getNewObject() override
double b
Definition: hdecay.h:120
float chiSquared() const
Definition: Trajectory.h:262
void getPeakOfLandau(TH1 *InputHisto, double *FitResults, double LowRange=0, double HighRange=5400)
SiStripCluster const & stereoCluster() const
edm::EventID id() const
Definition: EventBase.h:60
fixed size matrix
HLT enums.
GlobalVector globalMomentum() const
double a
Definition: hdecay.h:121
T get() const
Definition: EventSetup.h:68
const_iterator begin() const
first iterator over the map (read only)
bool isValid() const
Definition: ESHandle.h:45
long double T
T x() const
Definition: PV3DBase.h:62
unsigned int SubDet
Definition: APVGainStruct.h:12
T const * product() const
Definition: ESHandle.h:84
TimeValue_t value() const
Definition: Timestamp.h:56
const std::vector< uint8_t > & amplitudes() const
edm::Timestamp time() const
Definition: EventBase.h:61
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
const edm::Event * iEvent_
Definition: Run.h:44
unsigned int tecSide(const DetId &id) const
const SiStripApvGain::Range getRange(uint32_t detID) const
Definition: SiStripGain.h:70