#include <CalibCalorimetry/HcalStandardModules/interface/HFPreLightCal.h>
Public Member Functions | |
virtual void | analyze (const edm::Event &fEvent, const edm::EventSetup &fSetup) |
virtual void | beginJob (const edm::EventSetup &fSetup) |
virtual void | endJob (void) |
HFPreLightCal (const edm::ParameterSet &fConfiguration) | |
virtual | ~HFPreLightCal () |
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 |
Definition at line 11 of file HFPreLightCal.h.
HFPreLightCal::HFPreLightCal | ( | const edm::ParameterSet & | fConfiguration | ) |
Definition at line 38 of file HFPreLightCal.cc.
References edm::ParameterSet::getUntrackedParameter(), histfile, and textfile.
00038 { 00039 //std::string histfile = fConfiguration.getUntrackedParameter<string>("rootFile"); 00040 histfile = fConfiguration.getUntrackedParameter<string>("rootPreFile"); 00041 textfile = fConfiguration.getUntrackedParameter<string>("textPreFile"); 00042 }
HFPreLightCal::~HFPreLightCal | ( | ) | [virtual] |
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, GenMuonPlsPt100GeV_cfg::cout, HcalDetId::depth(), detId, lat::endl(), edm::EventID::event(), event_NN, edm::Event::getByType(), HcalForward, hts, edm::Event::id(), HFDataFrame::id(), HcalDetId::ieta(), HcalDetId::iphi(), edm::EventID::run(), run_NN, HFDataFrame::size(), and HcalDetId::subdet().
00124 { 00125 00126 // event ID 00127 edm::EventID eventId = fEvent.id(); 00128 int runNumber = eventId.run (); 00129 int eventNumber = eventId.event (); 00130 if (run_NN==0) run_NN=runNumber; 00131 event_NN++; 00132 if (verbose) std::cout << "========================================="<<std::endl 00133 << "run/event: "<<runNumber<<'/'<<eventNumber<<std::endl; 00134 00135 // HF PIN-diodes 00136 edm::Handle<HcalCalibDigiCollection> calib; 00137 fEvent.getByType(calib); 00138 if (verbose) std::cout<<"Analysis-> total CAL digis= "<<calib->size()<<std::endl; 00139 /* COMMENTED OUT by J. Mans (7-28-2008) as major changes needed with new Calib DetId 00140 00141 for (unsigned j = 0; j < calib->size (); ++j) { 00142 const HcalCalibDataFrame digi = (*calib)[j]; 00143 HcalElectronicsId elecId = digi.elecId(); 00144 HcalCalibDetId calibId = digi.id(); 00145 if (verbose) std::cout<<calibId.sectorString().c_str()<<" "<<calibId.rbx()<<" "<<elecId.fiberChanId()<<std::endl; 00146 int isector = calibId.rbx()-1; 00147 int ipin = elecId.fiberChanId(); 00148 int iside = -1; 00149 if (calibId.sectorString() == "HFP") iside = 0; 00150 else if (calibId.sectorString() == "HFM") iside = 4; 00151 00152 if (iside != -1) { 00153 for (int isample = 0; isample < digi.size(); ++isample) { 00154 int adc = digi[isample].adc(); 00155 int capid = digi[isample].capid (); 00156 double linear_ADC = digi[isample].nominal_fC(); 00157 if (verbose) std::cout<<"PIN linear_ADC = "<<linear_ADC<<std::endl; 00158 htspin[isector+iside][ipin]->Fill(isample,linear_ADC); 00159 } 00160 } 00161 } 00162 */ 00163 // HF 00164 edm::Handle<HFDigiCollection> hf_digi; 00165 fEvent.getByType(hf_digi); 00166 if (verbose) std::cout<<"Analysis-> total HF digis= "<<hf_digi->size()<<std::endl; 00167 00168 for (unsigned ihit = 0; ihit < hf_digi->size (); ++ihit) { 00169 const HFDataFrame& frame = (*hf_digi)[ihit]; 00170 HcalDetId detId = frame.id(); 00171 int ieta = detId.ieta(); 00172 int iphi = detId.iphi(); 00173 int depth = detId.depth(); 00174 if (verbose) std::cout <<"HF digi # " <<ihit<<": eta/phi/depth: " 00175 <<ieta<<'/'<<iphi<<'/'<< depth << std::endl; 00176 00177 if (ieta>0) ieta = ieta-29; 00178 else ieta = 13-ieta-29; 00179 00180 for (int isample = 0; isample < frame.size(); ++isample) { 00181 int adc = frame[isample].adc(); 00182 int capid = frame[isample].capid (); 00183 double linear_ADC = frame[isample].nominal_fC(); 00184 double nominal_fC = detId.subdet () == HcalForward ? 2.6 * linear_ADC : linear_ADC; 00185 00186 if (verbose) std::cout << "Analysis-> HF sample # " << isample 00187 << ", capid=" << capid 00188 << ": ADC=" << adc 00189 << ", linearized ADC=" << linear_ADC 00190 << ", nominal fC=" << nominal_fC << std::endl; 00191 00192 hts[ieta][(iphi-1)/2][depth-1]->Fill(isample,linear_ADC); 00193 } 00194 } 00195 return; 00196 }
void HFPreLightCal::beginJob | ( | const edm::EventSetup & | fSetup | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 48 of file HFPreLightCal.cc.
References GenMuonPlsPt100GeV_cfg::cout, lat::endl(), cmsRelvalreport::exit, histfile, hts, htsmax, htspin, htspinmax, i, j, k, mFile, NULL, textfile, and tFile.
00048 { 00049 00050 char htit[64]; 00051 std::cout<<std::endl<<"HFPreLightCal beginJob: --> "; 00052 00053 mFile = new TFile (histfile.c_str(),"RECREATE"); 00054 if ((tFile = fopen(textfile.c_str(),"w"))==NULL) { 00055 printf("\nNo textfile open\n\n"); 00056 std::cout<<"Problem with output Pre-textFILE => exit"<<std::endl; 00057 exit(1); 00058 } 00059 // Histos 00060 htsmax = new TH1F("htsmax","Max TS",100,0,10); 00061 htspinmax = new TH1F("htspinmax","Max TS PIN",100,0,10); 00062 // Channel-by-channel histos 00063 for (int i=0;i<13;i++) for (int j=0;j<36;j++) for (int k=0;k<2;k++) { 00064 if (i>10 && j%2==0) continue; 00065 sprintf(htit,"tspre_+%d_%d_%d",i+29,j*2+1,k+1); 00066 hts[i][j][k] = new TH1F(htit,htit,10,-0.5,9.5); // TimeSlices (pulse shape) 00067 sprintf(htit,"tspre_-%d_%d_%d",i+29,j*2+1,k+1); 00068 hts[i+13][j][k] = new TH1F(htit,htit,10,-0.5,9.5); // TimeSlices (pulse shape) 00069 } 00070 // PIN-diodes histos 00071 for (int i=0;i<4;i++) for (int j=0;j<3;j++) { 00072 sprintf(htit,"tspre_PIN%d_+Q%d",j+1,i+1); 00073 htspin[i][j] = new TH1F(htit,htit,10,-0.5,9.5); 00074 sprintf(htit,"tspre_PIN%d_-Q%d",j+1,i+1); 00075 htspin[i+4][j] = new TH1F(htit,htit,10,-0.5,9.5); 00076 } 00077 std::cout<<"histfile="<<histfile.c_str()<<" textfile="<<textfile.c_str()<<std::endl; 00078 return; 00079 }
Reimplemented from edm::EDAnalyzer.
Definition at line 81 of file HFPreLightCal.cc.
References cont, GenMuonPlsPt100GeV_cfg::cout, lat::endl(), event_NN, hts, htsmax, htspin, htspinmax, i, j, k, mFile, sum(), and tFile.
00082 { 00083 Double_t sum,cont; 00084 Int_t tsmax; 00085 00086 std::cout<<std::endl<<"HFPreLightCal endJob --> "; 00087 00088 for (int i=0;i<26;i++) for (int j=0;j<36;j++) for (int k=0;k<2;k++) { 00089 if (i>10 && i<13 && j%2==0) continue; 00090 if (i>23 && j%2==0) continue; 00091 sum=tsmax=0; 00092 for (int ii=1; ii<=10; ii++) { 00093 cont = hts[i][j][k]->GetBinContent(ii); 00094 if (ii<3) cont=cont-(hts[i][j][k]->GetBinContent(ii+4)+hts[i][j][k]->GetBinContent(ii+8))/2; 00095 else if (ii<5) cont=cont-hts[i][j][k]->GetBinContent(ii+4); 00096 else if (ii<7) cont=cont-(hts[i][j][k]->GetBinContent(ii-4)+hts[i][j][k]->GetBinContent(ii+4))/2; 00097 else if (ii<9) cont=cont-hts[i][j][k]->GetBinContent(ii-4); 00098 else cont=cont-(hts[i][j][k]->GetBinContent(ii-4)+hts[i][j][k]->GetBinContent(ii-8))/2; 00099 if (cont>sum) { 00100 sum = cont; 00101 tsmax=ii; 00102 } 00103 } 00104 htsmax->Fill(tsmax); 00105 if (i<13) fprintf(tFile," %d %d %d %d\n",i+29,j*2+1,k+1,tsmax); 00106 else fprintf(tFile," %d %d %d %d\n",13-i-29,j*2+1,k+1,tsmax); 00107 } 00108 00109 for (int i=0;i<8;i++) for (int j=0;j<3;j++) { 00110 sum=tsmax=0; 00111 tsmax = htspin[i][j]->GetMaximumBin(); 00112 htspinmax->Fill(tsmax); 00113 if (i<4) fprintf(tFile,"%d %d %d\n",i+1,j+1,tsmax); 00114 else fprintf(tFile,"%d %d %d\n",-i+3,j+1,tsmax); 00115 } 00116 00117 mFile->Write(); 00118 mFile->Close(); 00119 fclose(tFile); 00120 std::cout<<" Nevents = "<<event_NN<<std::endl; 00121 return; 00122 }
std::string HFPreLightCal::histfile [private] |
TH1F* HFPreLightCal::hts[26][36][2] [private] |
TH1F* HFPreLightCal::htsmax [private] |
TH1F* HFPreLightCal::htspin[8][3] [private] |
TH1F * HFPreLightCal::htspinmax [private] |
TFile* HFPreLightCal::mFile [private] |
std::string HFPreLightCal::textfile [private] |
FILE* HFPreLightCal::tFile [private] |