CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
HFPreLightCal Class Reference

#include <HFPreLightCal.h>

Inheritance diagram for HFPreLightCal:
edm::EDAnalyzer

Public Member Functions

virtual void analyze (const edm::Event &fEvent, const edm::EventSetup &fSetup)
 
virtual void beginJob ()
 
virtual void endJob (void)
 
 HFPreLightCal (const edm::ParameterSet &fConfiguration)
 
virtual ~HFPreLightCal ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Attributes

std::string histfile
 
TH1F * hts [26][36][2]
 
TH1F * htsmax
 
TH1F * htspin [8][3]
 
TH1F * htspinmax
 
TFile * mFile
 
std::string textfile
 
FILE * tFile
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Definition at line 11 of file HFPreLightCal.h.

Constructor & Destructor Documentation

HFPreLightCal::HFPreLightCal ( const edm::ParameterSet fConfiguration)

Definition at line 38 of file HFPreLightCal.cc.

References edm::ParameterSet::getUntrackedParameter().

38  {
39  //std::string histfile = fConfiguration.getUntrackedParameter<string>("rootFile");
40  histfile = fConfiguration.getUntrackedParameter<string>("rootPreFile");
41  textfile = fConfiguration.getUntrackedParameter<string>("textPreFile");
42 }
T getUntrackedParameter(std::string const &, T const &) const
std::string histfile
Definition: HFPreLightCal.h:26
std::string textfile
Definition: HFPreLightCal.h:27
HFPreLightCal::~HFPreLightCal ( )
virtual

Definition at line 44 of file HFPreLightCal.cc.

44  {
45  //delete mFile;
46 }

Member Function Documentation

void HFPreLightCal::analyze ( const edm::Event fEvent,
const edm::EventSetup fSetup 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 124 of file HFPreLightCal.cc.

References ecalMGPA::adc(), calib, gather_cfg::cout, HcalDetId::depth(), edm::EventID::event(), event_NN, HcalObjRepresent::Fill(), edm::Event::getByType(), HcalForward, HFDataFrame::id(), edm::EventBase::id(), HcalDetId::ieta(), HcalDetId::iphi(), edm::EventID::run(), run_NN, convertSQLiteXML::runNumber, HFDataFrame::size(), and HcalDetId::subdet().

124  {
125 
126  // event ID
127  edm::EventID eventId = fEvent.id();
128  int runNumber = eventId.run ();
129  int eventNumber = eventId.event ();
130  if (run_NN==0) run_NN=runNumber;
131  event_NN++;
132  if (verbose) std::cout << "========================================="<<std::endl
133  << "run/event: "<<runNumber<<'/'<<eventNumber<<std::endl;
134 
135  // HF PIN-diodes
137  fEvent.getByType(calib);
138  if (verbose) std::cout<<"Analysis-> total CAL digis= "<<calib->size()<<std::endl;
139  /* COMMENTED OUT by J. Mans (7-28-2008) as major changes needed with new Calib DetId
140 
141  for (unsigned j = 0; j < calib->size (); ++j) {
142  const HcalCalibDataFrame digi = (*calib)[j];
143  HcalElectronicsId elecId = digi.elecId();
144  HcalCalibDetId calibId = digi.id();
145  if (verbose) std::cout<<calibId.sectorString().c_str()<<" "<<calibId.rbx()<<" "<<elecId.fiberChanId()<<std::endl;
146  int isector = calibId.rbx()-1;
147  int ipin = elecId.fiberChanId();
148  int iside = -1;
149  if (calibId.sectorString() == "HFP") iside = 0;
150  else if (calibId.sectorString() == "HFM") iside = 4;
151 
152  if (iside != -1) {
153  for (int isample = 0; isample < digi.size(); ++isample) {
154  int adc = digi[isample].adc();
155  int capid = digi[isample].capid ();
156  double linear_ADC = digi[isample].nominal_fC();
157  if (verbose) std::cout<<"PIN linear_ADC = "<<linear_ADC<<std::endl;
158  htspin[isector+iside][ipin]->Fill(isample,linear_ADC);
159  }
160  }
161  }
162  */
163  // HF
165  fEvent.getByType(hf_digi);
166  if (verbose) std::cout<<"Analysis-> total HF digis= "<<hf_digi->size()<<std::endl;
167 
168  for (unsigned ihit = 0; ihit < hf_digi->size (); ++ihit) {
169  const HFDataFrame& frame = (*hf_digi)[ihit];
170  HcalDetId detId = frame.id();
171  int ieta = detId.ieta();
172  int iphi = detId.iphi();
173  int depth = detId.depth();
174  if (verbose) std::cout <<"HF digi # " <<ihit<<": eta/phi/depth: "
175  <<ieta<<'/'<<iphi<<'/'<< depth << std::endl;
176 
177  if (ieta>0) ieta = ieta-29;
178  else ieta = 13-ieta-29;
179 
180  for (int isample = 0; isample < frame.size(); ++isample) {
181  int adc = frame[isample].adc();
182  int capid = frame[isample].capid ();
183  double linear_ADC = frame[isample].nominal_fC();
184  double nominal_fC = detId.subdet () == HcalForward ? 2.6 * linear_ADC : linear_ADC;
185 
186  if (verbose) std::cout << "Analysis-> HF sample # " << isample
187  << ", capid=" << capid
188  << ": ADC=" << adc
189  << ", linearized ADC=" << linear_ADC
190  << ", nominal fC=" << nominal_fC << std::endl;
191 
192  hts[ieta][(iphi-1)/2][depth-1]->Fill(isample,linear_ADC);
193  }
194  }
195  return;
196 }
int adc(sample_type sample)
get the ADC sample (12 bits)
RunNumber_t run() const
Definition: EventID.h:42
EventNumber_t event() const
Definition: EventID.h:44
Int_t run_NN
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:32
Int_t event_NN
bool getByType(Handle< PROD > &result) const
Definition: Event.h:398
TH1F * hts[26][36][2]
Definition: HFPreLightCal.h:30
int depth() const
get the tower depth
Definition: HcalDetId.h:42
MVATrainerComputer * calib
Definition: MVATrainer.cc:64
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
int ieta() const
get the cell ieta
Definition: HcalDetId.h:38
int iphi() const
get the cell iphi
Definition: HcalDetId.h:40
int size() const
total number of samples in the digi
Definition: HFDataFrame.h:26
edm::EventID id() const
Definition: EventBase.h:56
tuple cout
Definition: gather_cfg.py:121
const HcalDetId & id() const
Definition: HFDataFrame.h:22
void HFPreLightCal::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 48 of file HFPreLightCal.cc.

References gather_cfg::cout, cmsRelvalreport::exit, i, j, gen::k, and NULL.

48  {
49 
50  char htit[64];
51  std::cout<<std::endl<<"HFPreLightCal beginJob: --> ";
52 
53  mFile = new TFile (histfile.c_str(),"RECREATE");
54  if ((tFile = fopen(textfile.c_str(),"w"))==NULL) {
55  printf("\nNo textfile open\n\n");
56  std::cout<<"Problem with output Pre-textFILE => exit"<<std::endl;
57  exit(1);
58  }
59  // Histos
60  htsmax = new TH1F("htsmax","Max TS",100,0,10);
61  htspinmax = new TH1F("htspinmax","Max TS PIN",100,0,10);
62  // Channel-by-channel histos
63  for (int i=0;i<13;i++) for (int j=0;j<36;j++) for (int k=0;k<2;k++) {
64  if (i>10 && j%2==0) continue;
65  sprintf(htit,"tspre_+%d_%d_%d",i+29,j*2+1,k+1);
66  hts[i][j][k] = new TH1F(htit,htit,10,-0.5,9.5); // TimeSlices (pulse shape)
67  sprintf(htit,"tspre_-%d_%d_%d",i+29,j*2+1,k+1);
68  hts[i+13][j][k] = new TH1F(htit,htit,10,-0.5,9.5); // TimeSlices (pulse shape)
69  }
70  // PIN-diodes histos
71  for (int i=0;i<4;i++) for (int j=0;j<3;j++) {
72  sprintf(htit,"tspre_PIN%d_+Q%d",j+1,i+1);
73  htspin[i][j] = new TH1F(htit,htit,10,-0.5,9.5);
74  sprintf(htit,"tspre_PIN%d_-Q%d",j+1,i+1);
75  htspin[i+4][j] = new TH1F(htit,htit,10,-0.5,9.5);
76  }
77  std::cout<<"histfile="<<histfile.c_str()<<" textfile="<<textfile.c_str()<<std::endl;
78  return;
79 }
int i
Definition: DBlmapReader.cc:9
std::string histfile
Definition: HFPreLightCal.h:26
#define NULL
Definition: scimark2.h:8
TH1F * hts[26][36][2]
Definition: HFPreLightCal.h:30
int j
Definition: DBlmapReader.cc:9
TH1F * htspin[8][3]
Definition: HFPreLightCal.h:32
int k[5][pyjets_maxn]
std::string textfile
Definition: HFPreLightCal.h:27
tuple cout
Definition: gather_cfg.py:121
TH1F * htspinmax
Definition: HFPreLightCal.h:31
void HFPreLightCal::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 81 of file HFPreLightCal.cc.

References cont, gather_cfg::cout, event_NN, i, j, and gen::k.

82 {
83  Double_t sum,cont;
84  Int_t tsmax;
85 
86  std::cout<<std::endl<<"HFPreLightCal endJob --> ";
87 
88  for (int i=0;i<26;i++) for (int j=0;j<36;j++) for (int k=0;k<2;k++) {
89  if (i>10 && i<13 && j%2==0) continue;
90  if (i>23 && j%2==0) continue;
91  sum=tsmax=0;
92  for (int ii=1; ii<=10; ii++) {
93  cont = hts[i][j][k]->GetBinContent(ii);
94  if (ii<3) cont=cont-(hts[i][j][k]->GetBinContent(ii+4)+hts[i][j][k]->GetBinContent(ii+8))/2;
95  else if (ii<5) cont=cont-hts[i][j][k]->GetBinContent(ii+4);
96  else if (ii<7) cont=cont-(hts[i][j][k]->GetBinContent(ii-4)+hts[i][j][k]->GetBinContent(ii+4))/2;
97  else if (ii<9) cont=cont-hts[i][j][k]->GetBinContent(ii-4);
98  else cont=cont-(hts[i][j][k]->GetBinContent(ii-4)+hts[i][j][k]->GetBinContent(ii-8))/2;
99  if (cont>sum) {
100  sum = cont;
101  tsmax=ii;
102  }
103  }
104  htsmax->Fill(tsmax);
105  if (i<13) fprintf(tFile," %d %d %d %d\n",i+29,j*2+1,k+1,tsmax);
106  else fprintf(tFile," %d %d %d %d\n",13-i-29,j*2+1,k+1,tsmax);
107  }
108 
109  for (int i=0;i<8;i++) for (int j=0;j<3;j++) {
110  sum=tsmax=0;
111  tsmax = htspin[i][j]->GetMaximumBin();
112  htspinmax->Fill(tsmax);
113  if (i<4) fprintf(tFile,"%d %d %d\n",i+1,j+1,tsmax);
114  else fprintf(tFile,"%d %d %d\n",-i+3,j+1,tsmax);
115  }
116 
117  mFile->Write();
118  mFile->Close();
119  fclose(tFile);
120  std::cout<<" Nevents = "<<event_NN<<std::endl;
121  return;
122 }
int i
Definition: DBlmapReader.cc:9
Int_t event_NN
TH1F * hts[26][36][2]
Definition: HFPreLightCal.h:30
int j
Definition: DBlmapReader.cc:9
TH1F * htspin[8][3]
Definition: HFPreLightCal.h:32
int k[5][pyjets_maxn]
int cont
tuple cout
Definition: gather_cfg.py:121
TH1F * htspinmax
Definition: HFPreLightCal.h:31

Member Data Documentation

std::string HFPreLightCal::histfile
private

Definition at line 26 of file HFPreLightCal.h.

TH1F* HFPreLightCal::hts[26][36][2]
private

Definition at line 30 of file HFPreLightCal.h.

TH1F* HFPreLightCal::htsmax
private

Definition at line 31 of file HFPreLightCal.h.

TH1F* HFPreLightCal::htspin[8][3]
private

Definition at line 32 of file HFPreLightCal.h.

TH1F * HFPreLightCal::htspinmax
private

Definition at line 31 of file HFPreLightCal.h.

TFile* HFPreLightCal::mFile
private

Definition at line 28 of file HFPreLightCal.h.

std::string HFPreLightCal::textfile
private

Definition at line 27 of file HFPreLightCal.h.

FILE* HFPreLightCal::tFile
private

Definition at line 29 of file HFPreLightCal.h.