00001
00002 #include "Validation/EcalRecHits/interface/EcalTBValidation.h"
00003 #include "DataFormats/EcalRecHit/interface/EcalUncalibratedRecHit.h"
00004 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00005 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
00006 #include "TBDataFormats/EcalTBObjects/interface/EcalTBHodoscopeRecInfo.h"
00007 #include "TBDataFormats/EcalTBObjects/interface/EcalTBTDCRecInfo.h"
00008 #include "TBDataFormats/EcalTBObjects/interface/EcalTBEventHeader.h"
00009 #include "DQMServices/Core/interface/DQMStore.h"
00010 #include "FWCore/ServiceRegistry/interface/Service.h"
00011
00012
00013 #include <iostream>
00014 #include <string>
00015
00016
00017 EcalTBValidation::EcalTBValidation( const edm::ParameterSet& config ) {
00018
00019 data_ = config.getUntrackedParameter<int>("data",-1000);
00020 xtalInBeam_ = config.getUntrackedParameter<int>("xtalInBeam",-1000);
00021
00022 digiCollection_ = config.getParameter<std::string>("digiCollection");
00023 digiProducer_ = config.getParameter<std::string>("digiProducer");
00024 hitCollection_ = config.getParameter<std::string>("hitCollection");
00025 hitProducer_ = config.getParameter<std::string>("hitProducer");
00026 hodoRecInfoCollection_ = config.getParameter<std::string>("hodoRecInfoCollection");
00027 hodoRecInfoProducer_ = config.getParameter<std::string>("hodoRecInfoProducer");
00028 tdcRecInfoCollection_ = config.getParameter<std::string>("tdcRecInfoCollection");
00029 tdcRecInfoProducer_ = config.getParameter<std::string>("tdcRecInfoProducer");
00030 eventHeaderCollection_ = config.getParameter<std::string>("eventHeaderCollection");
00031 eventHeaderProducer_ = config.getParameter<std::string>("eventHeaderProducer");
00032
00033
00034
00035
00036 verbose_ = config.getUntrackedParameter<bool>("verbose", false);
00037
00038 dbe_ = edm::Service<DQMStore>().operator->();
00039 if( dbe_ ) {
00040 if( verbose_ ) {
00041 dbe_->setVerbose(1);
00042 dbe_->showDirStructure();
00043 }
00044 else {
00045 dbe_->setVerbose(0);
00046 }
00047 }
00048
00049 meETBxib_ = 0;
00050 meETBampltdc_ = 0;
00051 meETBShape_ = 0;
00052 meETBhodoX_ = 0;
00053 meETBhodoY_ = 0;
00054 meETBe1x1_ = 0;
00055 meETBe3x3_ = 0;
00056 meETBe5x5_ = 0;
00057 meETBe1e9_ = 0;
00058 meETBe1e25_ = 0;
00059 meETBe9e25_ = 0;
00060 meETBe1x1_center_ = 0;
00061 meETBe3x3_center_ = 0;
00062 meETBe5x5_center_ = 0;
00063 meETBe1vsX_ = 0;
00064 meETBe1vsY_ = 0;
00065 meETBe1e9vsX_ = 0;
00066 meETBe1e9vsY_ = 0;
00067 meETBe1e25vsX_ = 0;
00068 meETBe1e25vsY_ = 0;
00069 meETBe9e25vsX_ = 0;
00070 meETBe9e25vsY_ = 0;
00071
00072 if( dbe_ ) {
00073
00074 std::string hname;
00075 dbe_->setCurrentFolder( "EcalRecHitsV/EcalTBValidationTask" );
00076
00077 hname = "xtal in beam position";
00078 meETBxib_ = dbe_->book2D( hname, hname, 85, 0., 85., 20,0., 20. );
00079 hname = "Max Amplitude vs TDC offset";
00080 meETBampltdc_ = dbe_->book2D( hname, hname, 100, 0., 1., 1000, 0., 4000. );
00081 hname = "Beam Profile X";
00082 meETBhodoX_ = dbe_->book1D( hname, hname, 100, -20., 20. );
00083 hname = "Beam Profile Y";
00084 meETBhodoY_ = dbe_->book1D( hname, hname, 100, -20., 20. );
00085 hname = "E1x1 energy";
00086 meETBe1x1_ = dbe_->book1D( hname, hname, 1000, 0., 4000. );
00087 hname = "E3x3 energy";
00088 meETBe3x3_ = dbe_->book1D( hname, hname, 1000, 0., 4000. );
00089 hname = "E5x5 energy";
00090 meETBe5x5_ = dbe_->book1D( hname, hname, 1000, 0., 4000. );
00091 hname = "E1x1 energy center";
00092 meETBe1x1_center_ = dbe_->book1D( hname, hname, 1000, 0., 4000. );
00093 hname = "E3x3 energy center";
00094 meETBe3x3_center_ = dbe_->book1D( hname, hname, 1000, 0., 4000. );
00095 hname = "E5x5 energy center";
00096 meETBe5x5_center_ = dbe_->book1D( hname, hname, 1000, 0., 4000. );
00097 hname = "E1 over E9 ratio";
00098 meETBe1e9_ = dbe_->book1D( hname, hname, 600, 0., 1.2 );
00099 hname = "E1 over E25 ratio";
00100 meETBe1e25_ = dbe_->book1D( hname, hname, 600, 0., 1.2 );
00101 hname = "E9 over E25 ratio";
00102 meETBe9e25_ = dbe_->book1D( hname, hname, 600, 0., 1.2 );
00103 hname = "E1 vs X";
00104 meETBe1vsX_ = dbe_->book2D( hname, hname, 80, -20, 20, 1000, 0., 4000. );
00105 hname = "E1 vs Y";
00106 meETBe1vsY_ = dbe_->book2D( hname, hname, 80, -20, 20, 1000, 0., 4000. );
00107 hname = "E1 over E9 vs X";
00108 meETBe1e9vsX_ = dbe_->book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
00109 hname = "E1 over E9 vs Y";
00110 meETBe1e9vsY_ = dbe_->book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
00111 hname = "E1 over E25 vs X";
00112 meETBe1e25vsX_ = dbe_->book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
00113 hname = "E1 over E25 vs Y";
00114 meETBe1e25vsY_ = dbe_->book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
00115 hname = "E9 over E25 vs X";
00116 meETBe9e25vsX_ = dbe_->book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
00117 hname = "E9 over E25 vs Y";
00118 meETBe9e25vsY_ = dbe_->book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
00119 hname = "Xtal in Beam Shape";
00120 meETBShape_ = dbe_->book2D( hname, hname, 250, 0, 10, 350, 0, 3500 );
00121 }
00122
00123 }
00124
00125
00126 EcalTBValidation::~EcalTBValidation(){}
00127
00128 void EcalTBValidation::beginJob() {}
00129
00130 void EcalTBValidation::endJob() {}
00131
00132 void EcalTBValidation::analyze( const edm::Event& event, const edm::EventSetup& setup ) {
00133
00134 using namespace edm;
00135 using namespace cms;
00136
00137
00138 const EBDigiCollection* theDigis=0;
00139 Handle<EBDigiCollection> pdigis;
00140 event.getByLabel(digiProducer_, digiCollection_, pdigis);
00141 if(pdigis.isValid()){
00142 theDigis = pdigis.product();
00143 }
00144 else {
00145 std::cerr << "Error! can't get the product " << digiCollection_.c_str() << std::endl;
00146 return;
00147 }
00148
00149
00150 const EBUncalibratedRecHitCollection* theHits=0;
00151 Handle<EBUncalibratedRecHitCollection> phits;
00152 event.getByLabel(hitProducer_, hitCollection_, phits);
00153 if(phits.isValid()){
00154 theHits = phits.product();
00155 }
00156 else {
00157 std::cerr << "Error! can't get the product " << hitCollection_.c_str() << std::endl;
00158 return;
00159 }
00160
00161
00162 const EcalTBHodoscopeRecInfo* theHodo=0;
00163 Handle<EcalTBHodoscopeRecInfo> pHodo;
00164 event.getByLabel(hodoRecInfoProducer_, hodoRecInfoCollection_, pHodo);
00165 if(pHodo.isValid()){
00166 theHodo = pHodo.product();
00167 }
00168 else{
00169 std::cerr << "Error! can't get the product " << hodoRecInfoCollection_.c_str() << std::endl;
00170 return;
00171 }
00172
00173
00174 const EcalTBTDCRecInfo* theTDC=0;
00175 Handle<EcalTBTDCRecInfo> pTDC;
00176 event.getByLabel(tdcRecInfoProducer_, tdcRecInfoCollection_, pTDC);
00177 if(pTDC.isValid()){
00178 theTDC = pTDC.product();
00179 }
00180 else{
00181 std::cerr << "Error! can't get the product " << tdcRecInfoCollection_.c_str() << std::endl;
00182 return;
00183 }
00184
00185
00186 const EcalTBEventHeader* evtHeader=0;
00187 Handle<EcalTBEventHeader> pEventHeader;
00188 event.getByLabel(eventHeaderProducer_ , pEventHeader);
00189 if(pEventHeader.isValid()){
00190 evtHeader = pEventHeader.product();
00191 }
00192 else{
00193 std::cerr << "Error! can't get the product " << eventHeaderProducer_.c_str() << std::endl;
00194 return;
00195 }
00196
00197
00198
00199
00200 EBDetId xtalInBeamId(1,xtalInBeam_, EBDetId::SMCRYSTALMODE);
00201 if (xtalInBeamId==EBDetId(0)){ return; }
00202 int xibEta = xtalInBeamId.ieta();
00203 int xibPhi = xtalInBeamId.iphi();
00204
00205
00206 if (data_ && (evtHeader->tableIsMoving())) return;
00207
00208
00209 EBDetId Xtals5x5[25];
00210 for (unsigned int icry=0;icry<25;icry++){
00211 unsigned int row = icry/5;
00212 unsigned int column = icry%5;
00213 int ieta = xtalInBeamId.ieta()+column-2;
00214 int iphi = xtalInBeamId.iphi()+row-2;
00215 if(EBDetId::validDetId(ieta, iphi)){
00216 EBDetId tempId(ieta, iphi,EBDetId::ETAPHIMODE);
00217 if (tempId.ism()==1)
00218 Xtals5x5[icry] = tempId;
00219 else
00220 Xtals5x5[icry] = EBDetId(0);
00221 } else {
00222 Xtals5x5[icry] = EBDetId(0);
00223 }
00224 }
00225
00226
00227 double ampl1x1 = 0.;
00228 double ampl3x3 = 0.;
00229 double ampl5x5 = 0.;
00230 for (unsigned int icry=0;icry<25;icry++) {
00231 if (!Xtals5x5[icry].null()){
00232 double theAmpl = (theHits->find(Xtals5x5[icry]))->amplitude();
00233 ampl5x5 += theAmpl;
00234 if (icry==12){ampl1x1 = theAmpl;}
00235 if (icry==6 || icry==7 || icry==8 || icry==11 || icry==12 || icry==13 || icry==16 || icry==17 || icry==18){ampl3x3 += theAmpl;}
00236 }}
00237
00238
00239
00240 double sampleSave[10];
00241 for(int ii=0; ii < 10; ++ii){ sampleSave[ii] = 0.0; }
00242 EBDigiCollection::const_iterator thisDigi = theDigis->find(xtalInBeamId);
00243
00244 double eMax = 0.;
00245 if (thisDigi != theDigis->end()){
00246 EBDataFrame myDigi = (*thisDigi);
00247 for (int sample=0; sample < myDigi.size(); ++sample){
00248 double analogSample = myDigi.sample(sample).adc();
00249 sampleSave[sample] = analogSample;
00250 if ( eMax < analogSample ) {
00251 eMax = analogSample;
00252
00253 }
00254 }
00255 }
00256
00257
00258 double xBeam = theHodo->posX();
00259 double yBeam = theHodo->posY();
00260
00261
00262 meETBxib_ -> Fill(xibEta, xibPhi);
00263 meETBhodoX_ -> Fill(xBeam);
00264 meETBhodoY_ -> Fill(yBeam);
00265 meETBampltdc_ -> Fill(theTDC->offset(),ampl1x1);
00266 meETBe1x1_ -> Fill(ampl1x1);
00267 meETBe3x3_ -> Fill(ampl3x3);
00268 meETBe5x5_ -> Fill(ampl5x5);
00269 meETBe1e9_ -> Fill(ampl1x1/ampl3x3);
00270 meETBe1e25_ -> Fill(ampl1x1/ampl5x5);
00271 meETBe9e25_ -> Fill(ampl3x3/ampl5x5);
00272 meETBe1vsX_ -> Fill(xBeam,ampl1x1);
00273 meETBe1vsY_ -> Fill(yBeam,ampl1x1);
00274 meETBe1e9vsX_ -> Fill(xBeam,ampl1x1/ampl3x3);
00275 meETBe1e9vsY_ -> Fill(yBeam,ampl1x1/ampl3x3);
00276 meETBe1e25vsX_ -> Fill(xBeam,ampl1x1/ampl5x5);
00277 meETBe1e25vsY_ -> Fill(yBeam,ampl1x1/ampl5x5);
00278 meETBe9e25vsX_ -> Fill(xBeam,ampl3x3/ampl5x5);
00279 meETBe9e25vsY_ -> Fill(yBeam,ampl3x3/ampl5x5);
00280
00281 for(int ii=0; ii < 10; ++ii){ meETBShape_->Fill(double(ii)+theTDC->offset(),sampleSave[ii]); }
00282
00283 if ( (fabs(xBeam)<2.5) && (fabs(yBeam)<2.5) ){
00284 meETBe1x1_center_ -> Fill(ampl1x1);
00285 meETBe3x3_center_ -> Fill(ampl3x3);
00286 meETBe5x5_center_ -> Fill(ampl5x5);
00287 }
00288 }
00289
00290