00001
00008
00009
00010
00011
00012
00013 #include "RecoTBCalo/EcalSimpleTBAnalysis/interface/EcalSimple2007H4TBAnalyzer.h"
00014 #include "DataFormats/EcalRecHit/interface/EcalUncalibratedRecHit.h"
00015 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00016 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
00017 #include "TBDataFormats/EcalTBObjects/interface/EcalTBHodoscopeRecInfo.h"
00018 #include "TBDataFormats/EcalTBObjects/interface/EcalTBTDCRecInfo.h"
00019 #include "TBDataFormats/EcalTBObjects/interface/EcalTBEventHeader.h"
00020
00021 #include "FWCore/Framework/interface/ESHandle.h"
00022
00023 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00024 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00025 #include "Geometry/CaloGeometry/interface/TruncatedPyramid.h"
00026
00027
00028
00029 #include "TFile.h"
00030 #include "TH1.h"
00031 #include "TH2.h"
00032 #include "TF1.h"
00033
00034 #include <iostream>
00035 #include <string>
00036 #include <stdexcept>
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051 EcalSimple2007H4TBAnalyzer::EcalSimple2007H4TBAnalyzer( const edm::ParameterSet& iConfig ) : xtalInBeam_(0)
00052
00053 {
00054
00055 rootfile_ = iConfig.getUntrackedParameter<std::string>("rootfile","ecalSimpleTBanalysis.root");
00056 digiCollection_ = iConfig.getParameter<std::string>("digiCollection");
00057 digiProducer_ = iConfig.getParameter<std::string>("digiProducer");
00058 hitCollection_ = iConfig.getParameter<std::string>("hitCollection");
00059 hitProducer_ = iConfig.getParameter<std::string>("hitProducer");
00060 hodoRecInfoCollection_ = iConfig.getParameter<std::string>("hodoRecInfoCollection");
00061 hodoRecInfoProducer_ = iConfig.getParameter<std::string>("hodoRecInfoProducer");
00062 tdcRecInfoCollection_ = iConfig.getParameter<std::string>("tdcRecInfoCollection");
00063 tdcRecInfoProducer_ = iConfig.getParameter<std::string>("tdcRecInfoProducer");
00064 eventHeaderCollection_ = iConfig.getParameter<std::string>("eventHeaderCollection");
00065 eventHeaderProducer_ = iConfig.getParameter<std::string>("eventHeaderProducer");
00066
00067
00068 std::cout << "EcalSimple2007H4TBAnalyzer: fetching hitCollection: " << hitCollection_.c_str()
00069 << " produced by " << hitProducer_.c_str() << std::endl;
00070
00071 }
00072
00073
00074
00075 EcalSimple2007H4TBAnalyzer::~EcalSimple2007H4TBAnalyzer()
00076
00077 {
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101 }
00102
00103
00104 void
00105 EcalSimple2007H4TBAnalyzer::beginRun(edm::Run const &, edm::EventSetup const& iSetup) {
00106
00107
00108 edm::ESHandle<CaloGeometry> pG;
00109 iSetup.get<CaloGeometryRecord>().get(pG);
00110
00111
00112 theTBGeometry_ = &(*pG);
00113
00114
00115
00116
00117
00118
00119 h_ampltdc = new TH2F("h_ampltdc","Max Amplitude vs TDC offset", 100,0.,1.,1000, 0., 4000.);
00120
00121
00122 h_tableIsMoving = new TH1F("h_tableIsMoving","TableIsMoving", 100000, 0., 100000.);
00123
00124 h_e1x1 = new TH1F("h_e1x1","E1x1 energy", 1000, 0., 4000.);
00125 h_e3x3 = new TH1F("h_e3x3","E3x3 energy", 1000, 0., 4000.);
00126 h_e5x5 = new TH1F("h_e5x5","E5x5 energy", 1000, 0., 4000.);
00127
00128 h_e1x1_center = new TH1F("h_e1x1_center","E1x1 energy", 1000, 0., 4000.);
00129 h_e3x3_center = new TH1F("h_e3x3_center","E3x3 energy", 1000, 0., 4000.);
00130 h_e5x5_center = new TH1F("h_e5x5_center","E5x5 energy", 1000, 0., 4000.);
00131
00132 h_e1e9 = new TH1F("h_e1e9","E1/E9 ratio", 600, 0., 1.2);
00133 h_e1e25 = new TH1F("h_e1e25","E1/E25 ratio", 600, 0., 1.2);
00134 h_e9e25 = new TH1F("h_e9e25","E9/E25 ratio", 600, 0., 1.2);
00135
00136 h_S6 = new TH1F("h_S6","Amplitude S6", 1000, 0., 4000.);
00137
00138 h_bprofx = new TH1F("h_bprofx","Beam Profile X",100,-20.,20.);
00139 h_bprofy = new TH1F("h_bprofy","Beam Profile Y",100,-20.,20.);
00140
00141 h_qualx = new TH1F("h_qualx","Beam Quality X",5000,0.,5.);
00142 h_qualy = new TH1F("h_qualy","Beam Quality X",5000,0.,5.);
00143
00144 h_slopex = new TH1F("h_slopex","Beam Slope X",500, -5e-4 , 5e-4 );
00145 h_slopey = new TH1F("h_slopey","Beam Slope Y",500, -5e-4 , 5e-4 );
00146
00147 char hname[50];
00148 char htitle[50];
00149 for (unsigned int icry=0;icry<25;icry++)
00150 {
00151 sprintf(hname,"h_mapx_%d",icry);
00152 sprintf(htitle,"Max Amplitude vs X %d",icry);
00153 h_mapx[icry] = new TH2F(hname,htitle,80,-20,20,1000,0.,4000.);
00154 sprintf(hname,"h_mapy_%d",icry);
00155 sprintf(htitle,"Max Amplitude vs Y %d",icry);
00156 h_mapy[icry] = new TH2F(hname,htitle,80,-20,20,1000,0.,4000.);
00157 }
00158
00159 h_e1e9_mapx = new TH2F("h_e1e9_mapx","E1/E9 vs X",80,-20,20,600,0.,1.2);
00160 h_e1e9_mapy = new TH2F("h_e1e9_mapy","E1/E9 vs Y",80,-20,20,600,0.,1.2);
00161
00162 h_e1e25_mapx = new TH2F("h_e1e25_mapx","E1/E25 vs X",80,-20,20,600,0.,1.2);
00163 h_e1e25_mapy = new TH2F("h_e1e25_mapy","E1/E25 vs Y",80,-20,20,600,0.,1.2);
00164
00165 h_e9e25_mapx = new TH2F("h_e9e25_mapx","E9/E25 vs X",80,-20,20,600,0.,1.2);
00166 h_e9e25_mapy = new TH2F("h_e9e25_mapy","E9/E25 vs Y",80,-20,20,600,0.,1.2);
00167
00168 h_Shape_ = new TH2F("h_Shape_","Xtal in Beam Shape",250,0,10,350,0,3500);
00169
00170 }
00171
00172
00173 void
00174 EcalSimple2007H4TBAnalyzer::endJob() {
00175
00176
00177 TFile f(rootfile_.c_str(),"RECREATE");
00178
00179
00180 h_ampltdc->Write();
00181
00182
00183 h_e1x1->Write();
00184 h_e3x3->Write();
00185 h_e5x5->Write();
00186
00187 h_e1x1_center->Write();
00188 h_e3x3_center->Write();
00189 h_e5x5_center->Write();
00190
00191 h_e1e9->Write();
00192 h_e1e25->Write();
00193 h_e9e25->Write();
00194
00195 h_S6->Write();
00196 h_bprofx->Write();
00197 h_bprofy->Write();
00198
00199 h_qualx->Write();
00200 h_qualy->Write();
00201
00202 h_slopex->Write();
00203 h_slopey->Write();
00204
00205 h_Shape_->Write();
00206
00207 for (unsigned int icry=0;icry<25;icry++)
00208 {
00209 h_mapx[icry]->Write();
00210 h_mapy[icry]->Write();
00211 }
00212
00213 h_e1e9_mapx->Write();
00214 h_e1e9_mapy->Write();
00215
00216 h_e1e25_mapx->Write();
00217 h_e1e25_mapy->Write();
00218
00219 h_e9e25_mapx->Write();
00220 h_e9e25_mapy->Write();
00221
00222 h_tableIsMoving->Write();
00223
00224 f.Close();
00225 }
00226
00227
00228
00229
00230
00231
00232 void
00233 EcalSimple2007H4TBAnalyzer::analyze( edm::Event const & iEvent, edm::EventSetup const & iSetup ) {
00234
00235
00236 using namespace edm;
00237 using namespace cms;
00238
00239
00240
00241 Handle<EEDigiCollection> pdigis;
00242 const EEDigiCollection* digis=0;
00243
00244 iEvent.getByLabel( digiProducer_, digiCollection_,pdigis);
00245 if ( pdigis.isValid() ) {
00246 digis = pdigis.product();
00247
00248 } else {
00249 edm::LogError("EcalSimple2007H4TBAnalyzerError") << "Error! can't get the product " << digiCollection_;
00250 }
00251
00252
00253 Handle<EEUncalibratedRecHitCollection> phits;
00254 const EEUncalibratedRecHitCollection* hits=0;
00255
00256 iEvent.getByLabel( hitProducer_, hitCollection_,phits);
00257 if (phits.isValid()) {
00258 hits = phits.product();
00259
00260 } else {
00261 edm::LogError("EcalSimple2007H4TBAnalyzerError") << "Error! can't get the product " << hitCollection_;
00262 }
00263
00264 Handle<EcalTBHodoscopeRecInfo> pHodo;
00265 const EcalTBHodoscopeRecInfo* recHodo=0;
00266
00267 iEvent.getByLabel( hodoRecInfoProducer_, hodoRecInfoCollection_, pHodo);
00268 if ( pHodo.isValid() ) {
00269 recHodo = pHodo.product();
00270 } else {
00271 edm::LogError("EcalSimple2007H4TBAnalyzerError") << "Error! can't get the product " << hodoRecInfoCollection_;
00272 }
00273
00274 Handle<EcalTBTDCRecInfo> pTDC;
00275 const EcalTBTDCRecInfo* recTDC=0;
00276
00277 iEvent.getByLabel( tdcRecInfoProducer_, tdcRecInfoCollection_, pTDC);
00278 if ( pTDC.isValid() ) {
00279 recTDC = pTDC.product();
00280 } else {
00281 edm::LogError("EcalSimple2007H4TBAnalyzerError") << "Error! can't get the product " << tdcRecInfoCollection_;
00282 }
00283
00284 Handle<EcalTBEventHeader> pEventHeader;
00285 const EcalTBEventHeader* evtHeader=0;
00286
00287 iEvent.getByLabel( eventHeaderProducer_ , pEventHeader );
00288 if ( pEventHeader.isValid() ) {
00289 evtHeader = pEventHeader.product();
00290 } else {
00291 edm::LogError("EcalSimple2007H4TBAnalyzerError") << "Error! can't get the product " << eventHeaderProducer_;
00292 }
00293
00294
00295 if (!hits)
00296 {
00297
00298 return;
00299 }
00300
00301 if (!recTDC)
00302 {
00303
00304 return;
00305 }
00306
00307 if (!recHodo)
00308 {
00309
00310 return;
00311 }
00312
00313 if (!evtHeader)
00314 {
00315
00316 return;
00317 }
00318
00319 if (hits->size() == 0)
00320 {
00321
00322 return;
00323 }
00324
00325
00326 if (evtHeader->tableIsMoving())
00327 h_tableIsMoving->Fill(evtHeader->eventNumber());
00328
00329
00330
00331
00332
00333
00334
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344 h_S6->Fill(evtHeader->S6ADC());
00345
00346 if (xtalInBeamTmp.null())
00347 {
00348 xtalInBeamTmp = EBDetId(1,evtHeader->crystalInBeam(),EBDetId::SMCRYSTALMODE);
00349 xtalInBeam_ = EEDetId( 35 - ((xtalInBeamTmp.ic()-1)%20) ,int(int(xtalInBeamTmp.ic())/int(20))+51, -1);
00350 std::cout<< "Xtal In Beam is " << xtalInBeam_.ic() << xtalInBeam_ << std::endl;
00351 for (unsigned int icry=0;icry<25;icry++)
00352 {
00353 unsigned int row = icry / 5;
00354 unsigned int column= icry %5;
00355 int ix=xtalInBeam_.ix()+row-2;
00356 int iy=xtalInBeam_.iy()+column-2;
00357 EEDetId tempId(ix, iy, xtalInBeam_.zside());
00358
00359 if (tempId.ix()<16 || tempId.ix()>35 || tempId.iy()<51 || tempId.iy()>75)
00360 Xtals5x5[icry]=EEDetId(0);
00361 else
00362 {
00363 Xtals5x5[icry]=tempId;
00364 const CaloCellGeometry* cell=theTBGeometry_->getGeometry(Xtals5x5[icry]);
00365 if (!cell)
00366 continue;
00367 const TruncatedPyramid* tp ( dynamic_cast<const TruncatedPyramid*>(cell) ) ;
00368 std::cout << "** Xtal in the matrix **** row " << row << ", column " << column << ", xtal " << Xtals5x5[icry] << " Position " << tp->getPosition(0.) << std::endl;
00369 }
00370 }
00371 }
00372 else
00373 if (xtalInBeamTmp != EBDetId(1,evtHeader->crystalInBeam(),EBDetId::SMCRYSTALMODE))
00374 return;
00375
00376
00377 if (evtHeader->tableIsMoving())
00378 {
00379 std::cout << "Table is moving" << std::endl;
00380 return;
00381 }
00382
00383
00384
00385
00386 EEDetId maxHitId(0);
00387 float maxHit= -999999.;
00388 for(EEUncalibratedRecHitCollection::const_iterator ithit = hits->begin(); ithit != hits->end(); ++ithit)
00389 {
00390 if (ithit->amplitude()>=maxHit)
00391 {
00392 maxHit=ithit->amplitude();
00393 maxHitId=ithit->id();
00394 }
00395
00396 }
00397 if (maxHitId==EEDetId(0))
00398 {
00399 std::cout << "No maxHit found" << std::endl;
00400 return;
00401 }
00402
00403
00404
00405 double samples_save[10]; for(int i=0; i < 10; ++i) samples_save[i]=0.0;
00406
00407 double eMax = 0.;
00408 for ( EEDigiCollection::const_iterator digiItr= digis->begin();digiItr != digis->end();
00409 ++digiItr )
00410 {
00411 if ( EEDetId((*digiItr).id()) != xtalInBeam_ )
00412 continue;
00413
00414 EEDataFrame myDigi = (*digiItr);
00415 for (int sample = 0; sample < myDigi.size(); ++sample)
00416 {
00417 double analogSample = myDigi.sample(sample).adc();
00418 samples_save[sample] = analogSample;
00419
00420 if ( eMax < analogSample )
00421 {
00422 eMax = analogSample;
00423 }
00424 }
00425
00426 }
00427
00428 for(int i =0; i < 10; ++i)
00429 h_Shape_->Fill(double(i)+recTDC->offset(),samples_save[i]);
00430
00431
00432
00433
00434 double amplitude[25];
00435 double amplitude3x3=0;
00436 double amplitude5x5=0;
00437 for (unsigned int icry=0;icry<25;icry++)
00438 {
00439 if (!Xtals5x5[icry].null())
00440 {
00441 amplitude[icry]=(hits->find(Xtals5x5[icry]))->amplitude();
00442 amplitude5x5 += amplitude[icry];
00443
00444 if ( icry == 6 || icry == 7 || icry == 8 ||
00445 icry == 11 || icry == 12 || icry ==13 ||
00446 icry == 16 || icry == 17 || icry ==18 )
00447 {
00448 amplitude3x3+=amplitude[icry];
00449 }
00450 }
00451 }
00452
00453
00454 h_e1x1->Fill(amplitude[12]);
00455 h_e3x3->Fill(amplitude3x3);
00456 h_e5x5->Fill(amplitude5x5);
00457
00458 h_e1e9->Fill(amplitude[12]/amplitude3x3);
00459 h_e1e25->Fill(amplitude[12]/amplitude5x5);
00460 h_e9e25->Fill(amplitude3x3/amplitude5x5);
00461
00462
00463 if (recTDC)
00464 h_ampltdc->Fill(recTDC->offset(),amplitude[12]);
00465
00466
00467 if (recHodo)
00468 {
00469 float x=recHodo->posX();
00470 float y=recHodo->posY();
00471 float xslope=recHodo->slopeX();
00472 float yslope=recHodo->slopeY();
00473 float xqual=recHodo->qualX();
00474 float yqual=recHodo->qualY();
00475
00476
00477 h_bprofx->Fill(x);
00478 h_bprofy->Fill(y);
00479 h_qualx->Fill(xqual);
00480 h_qualy->Fill(yqual);
00481 h_slopex->Fill(xslope);
00482 h_slopey->Fill(yslope);
00483
00484
00485
00486
00487 if ( fabs(x + 2.5) < 2.5 && fabs(y + 0.5) < 2.5)
00488 {
00489 h_e1x1_center->Fill(amplitude[12]);
00490 h_e3x3_center->Fill(amplitude3x3);
00491 h_e5x5_center->Fill(amplitude5x5);
00492 }
00493
00494 for (unsigned int icry=0;icry<25;icry++)
00495 {
00496 h_mapx[icry]->Fill(x,amplitude[icry]);
00497 h_mapy[icry]->Fill(y,amplitude[icry]);
00498 }
00499
00500 h_e1e9_mapx->Fill(x,amplitude[12]/amplitude3x3);
00501 h_e1e9_mapy->Fill(y,amplitude[12]/amplitude3x3);
00502
00503 h_e1e25_mapx->Fill(x,amplitude[12]/amplitude5x5);
00504 h_e1e25_mapy->Fill(y,amplitude[12]/amplitude5x5);
00505
00506 h_e9e25_mapx->Fill(x,amplitude3x3/amplitude5x5);
00507 h_e9e25_mapy->Fill(y,amplitude3x3/amplitude5x5);
00508 }
00509
00510 }
00511
00512