00001 #include "Calibration/EcalCalibAlgos/src/PhiSymmetryCalibration_step2_SM.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00004 #include "CondFormats/EcalObjects/interface/EcalIntercalibConstants.h"
00005 #include "CondTools/Ecal/interface/EcalIntercalibConstantsXMLTranslator.h"
00006 #include "CondFormats/DataRecord/interface/EcalIntercalibConstantsRcd.h"
00007 #include "CondFormats/DataRecord/interface/EcalChannelStatusRcd.h"
00008 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00009 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00010
00011 #include "TH2F.h"
00012
00013 #include "TH1F.h"
00014 #include "TFile.h"
00015
00016 #include <fstream>
00017 #include "boost/filesystem/operations.hpp"
00018
00019 using namespace std;
00020
00021
00022
00023 PhiSymmetryCalibration_step2_SM::~PhiSymmetryCalibration_step2_SM(){}
00024
00025
00026 PhiSymmetryCalibration_step2_SM::PhiSymmetryCalibration_step2_SM(const edm::ParameterSet& iConfig){
00027
00028 statusThreshold_ =
00029 iConfig.getUntrackedParameter<int>("statusThreshold",0);
00030 have_initial_miscalib_=
00031 iConfig.getUntrackedParameter<bool>("haveInitialMiscalib",false);
00032 initialmiscalibfile_=
00033 iConfig.getUntrackedParameter<std::string>("initialmiscalibfile",
00034 "InitialMiscalib.xml");
00035 oldcalibfile_=
00036 iConfig.getUntrackedParameter<std::string>("oldcalibfile",
00037 "EcalIntercalibConstants.xml");
00038 reiteration_ = iConfig.getUntrackedParameter<bool>("reiteration",false);
00039 firstpass_=true;
00040 }
00041
00042 void PhiSymmetryCalibration_step2_SM::analyze( const edm::Event& ev,
00043 const edm::EventSetup& se){
00044
00045 if (firstpass_) {
00046 setUp(se);
00047 firstpass_=false;
00048 }
00049 }
00050
00051 void PhiSymmetryCalibration_step2_SM::setUp(const edm::EventSetup& se){
00052
00053 edm::ESHandle<EcalChannelStatus> chStatus;
00054 se.get<EcalChannelStatusRcd>().get(chStatus);
00055
00056 edm::ESHandle<CaloGeometry> geoHandle;
00057 se.get<CaloGeometryRecord>().get(geoHandle);
00058
00059 barrelCells = geoHandle->getValidDetIds(DetId::Ecal, EcalBarrel);
00060 endcapCells = geoHandle->getValidDetIds(DetId::Ecal, EcalEndcap);
00061
00062 e_.setup(&(*geoHandle), &(*chStatus), statusThreshold_);
00063
00064
00065
00066 for (int sign=0; sign<kSides; sign++) {
00067 for (int ieta=0; ieta<kBarlRings; ieta++) {
00068 for (int iphi=0; iphi<kBarlWedges; iphi++) {
00069 int iphi_r=int(iphi/nscx);
00070 if( !e_.goodCell_barl[ieta][iphi][sign] ){
00071 nBads_barl_SM_[ieta][iphi_r][sign]++;
00072
00073 }
00074 }
00075 }
00076 }
00077
00078
00079
00081 if (have_initial_miscalib_){
00082
00083 EcalCondHeader h;
00084 namespace fs = boost::filesystem;
00085 fs::path p(initialmiscalibfile_.c_str(),fs::native);
00086 if (!fs::exists(p)) edm::LogError("PhiSym") << "File not found: "
00087 << initialmiscalibfile_ <<endl;
00088
00089 int ret=
00090 EcalIntercalibConstantsXMLTranslator::readXML(initialmiscalibfile_,h,miscalib_);
00091 if (ret) edm::LogError("PhiSym")<<"Error reading XML files"<<endl;
00092 } else {
00093
00094 for (vector<DetId>::iterator it=barrelCells.begin(); it!=barrelCells.end(); ++it){
00095 miscalib_[*it]=1;
00096 }
00097
00098 for (vector<DetId>::iterator it=endcapCells.begin(); it!=endcapCells.end(); ++it){
00099 miscalib_[*it]=1;
00100
00101 }
00102 }
00103
00104
00105
00106 if (reiteration_){
00107
00108
00109 EcalCondHeader h;
00110 namespace fs = boost::filesystem;
00111 fs::path p(oldcalibfile_.c_str(),fs::native);
00112 if (!fs::exists(p)) edm::LogError("PhiSym") << "File not found: "
00113 << oldcalibfile_ <<endl;
00114
00115 int ret=
00116 EcalIntercalibConstantsXMLTranslator::readXML(oldcalibfile_,h,
00117 oldCalibs_);
00118
00119 if (ret) edm::LogError("PhiSym")<<"Error reading XML files"<<endl;;
00120
00121 } else {
00122
00123 for (vector<DetId>::iterator it=barrelCells.begin();
00124 it!=barrelCells.end(); ++it)
00125 oldCalibs_[*it]=1;
00126
00127
00128 for (vector<DetId>::iterator it=endcapCells.begin();
00129 it!=endcapCells.end(); ++it)
00130 oldCalibs_[*it]=1;
00131
00132
00133 }
00134
00135 }
00136
00137
00138 void PhiSymmetryCalibration_step2_SM::beginJob(){
00139
00140
00141 for (int ieta=0; ieta<kBarlRings; ieta++) {
00142 for (int iphi=0; iphi<kBarlWedges; iphi++) {
00143 for (int sign=0; sign<kSides; sign++) {
00144 int iphi_r=int(iphi/nscx);
00145
00146 etsum_barl_[ieta][iphi][sign]=0.;
00147 nhits_barl_[ieta][iphi][sign]=0;
00148 esum_barl_[ieta][iphi][sign]=0.;
00149 etsum_barl_SM_[ieta][iphi_r][sign]=0;
00150 nBads_barl_SM_[ieta][iphi_r][sign]=0;
00151 epsilon_M_barl_SM_[ieta][iphi_r][sign]=0;
00152 }
00153 }
00154 etsumMean_barl_SM_[ieta]=0.;
00155 }
00156
00157 for (int ix=0; ix<kEndcWedgesX; ix++) {
00158 for (int iy=0; iy<kEndcWedgesY; iy++) {
00159 for (int sign=0; sign<kSides; sign++) {
00160 etsum_endc_[ix][iy][sign]=0.;
00161 nhits_endc_[ix][iy][sign]=0;
00162 esum_endc_[ix][iy][sign]=0.;
00163
00164 }
00165 }
00166 }
00167
00168 readEtSums();
00169 setupResidHistos();
00170 }
00171
00172 void PhiSymmetryCalibration_step2_SM::endJob(){
00173
00174 if (firstpass_) {
00175 edm::LogError("PhiSym")<< "Must process at least one event-Exiting" <<endl;
00176 return;
00177
00178 }
00179
00180
00181
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203 fillHistos();
00204
00205
00206 std::ofstream etsumMean_barl_out("etsumMean_barl.dat",ios::out);
00207 for (int ieta=0; ieta<kBarlRings; ieta++) {
00208 etsumMean_barl_out << ieta << " " << etsumMean_barl_[ieta] << endl;
00209 }
00210 etsumMean_barl_out.close();
00211
00212 std::ofstream etsumMean_endc_out("etsumMean_endc.dat",ios::out);
00213 for (int ring=0; ring<kEndcEtaRings; ring++) {
00214 etsumMean_endc_out << e_.cellPos_[ring][50].eta() << " " << etsumMean_endc_[ring] << endl;
00215 }
00216 etsumMean_endc_out.close();
00217
00218
00219
00220 for (int ieta=0; ieta<kBarlRings; ieta++) {
00221 for (int iphi=0; iphi<kBarlWedges; iphi++) {
00222 for (int sign=0; sign<kSides; sign++) {
00223
00224
00225
00226 int iphi_r = int(iphi/nscx);
00227
00228
00229
00230
00231
00232 float epsilon_T_SM =
00233 etsum_barl_SM_[ieta][iphi_r][sign] /etsumMean_barl_SM_[ieta] -1.;
00234
00235 epsilon_M_barl_SM_[ieta][iphi_r][sign] = epsilon_T_SM/k_barl_[ieta];
00236
00237 if(e_.goodCell_barl[ieta][iphi][sign]){
00238 float etsum = etsum_barl_[ieta][iphi][sign];
00239 float epsilon_T = (etsum/etsumMean_barl_[ieta]) - 1.;
00240 rawconst_barl[ieta][iphi][sign] = epsilon_T + 1.;
00241 epsilon_M_barl[ieta][iphi][sign] = epsilon_T/k_barl_[ieta];
00242
00243
00244 } else {
00245 rawconst_barl[ieta][iphi][sign] = 1.;
00246 epsilon_M_barl[ieta][iphi][sign] = 0.;
00247 }
00248 }
00249 }
00250 }
00251
00252
00253 for (int ix=0; ix<kEndcWedgesX; ix++) {
00254 for (int iy=0; iy<kEndcWedgesY; iy++) {
00255 for (int sign=0; sign<kSides; sign++) {
00256 int ring = e_.endcapRing_[ix][iy];
00257 if (ring!=-1 && e_.goodCell_endc[ix][iy][sign]) {
00258 float etsum = etsum_endc_[ix][iy][sign];
00259 float epsilon_T = (etsum/etsumMean_endc_[ring]) - 1.;
00260 rawconst_endc[ix][iy][sign] = epsilon_T + 1.;
00261 epsilon_M_endc[ix][iy][sign] = epsilon_T/k_endc_[ring];
00262 } else {
00263 epsilon_M_endc[ix][iy][0] = 0.;
00264 epsilon_M_endc[ix][iy][1] = 0.;
00265 rawconst_endc[ix][iy][0] = 1.;
00266 rawconst_endc[ix][iy][1] = 1.;
00267 }
00268 }
00269 }
00270 }
00271
00272
00273
00274
00275 fstream scfile("sccalibration.dat",ios::out);
00276 for (int ieta =0; ieta< kBarlRings; ++ieta){
00277 for (int iphi_r =0; iphi_r< int(kBarlWedges/nscx);++iphi_r){
00278 for (int sign=0; sign<kSides; sign++) {
00279 scfile<< ieta << " " << iphi_r << " " <<sign << " "
00280 << 1/(1+epsilon_M_barl_SM_[ieta][iphi_r][sign] )<< std::endl;
00281 }
00282 }
00283 }
00284
00285 std::string newcalibfile("EcalIntercalibConstants_new.xml");
00286
00287
00288
00289 TFile ehistof("ehistos.root","recreate");
00290
00291 TH1D ebhisto("eb","eb",100, 0.,2.);
00292
00293 std::vector<DetId>::const_iterator barrelIt=barrelCells.begin();
00294 for (; barrelIt!=barrelCells.end(); barrelIt++) {
00295 EBDetId eb(*barrelIt);
00296 int ieta = abs(eb.ieta())-1;
00297 int iphi = eb.iphi()-1;
00298 int sign = eb.zside()>0 ? 1 : 0;
00299
00302 newCalibs_[eb] = oldCalibs_[eb]/(1+epsilon_M_barl[ieta][iphi][sign]);
00303
00304 if(e_.goodCell_barl[ieta][iphi][sign]){
00305
00306 ebhisto.Fill(newCalibs_[eb]);
00307
00309 miscal_resid_barl_histos[ieta]->Fill(miscalib_[eb]*newCalibs_[eb]);
00310 correl_barl_histos[ieta]->Fill(miscalib_[eb],newCalibs_[eb]);
00311
00312 }
00313
00314 }
00315
00316 TH1D eehisto("ee","ee",100, 0.,2.);
00317 std::vector<DetId>::const_iterator endcapIt=endcapCells.begin();
00318
00319 for (; endcapIt!=endcapCells.end(); endcapIt++) {
00320 EEDetId ee(*endcapIt);
00321 int ix = ee.ix()-1;
00322 int iy = ee.iy()-1;
00323 int sign = ee.zside()>0 ? 1 : 0;
00324
00325 newCalibs_[ee] = oldCalibs_[ee]/(1+epsilon_M_endc[ix][iy][sign]);
00326
00327
00328 if(e_.goodCell_endc[ix][iy][sign]){
00329
00330 eehisto.Fill(newCalibs_[ee]);
00331 miscal_resid_endc_histos[e_.endcapRing_[ix][iy]]->Fill(miscalib_[ee]*
00332 newCalibs_[ee]);;
00333
00334 correl_endc_histos[e_.endcapRing_[ix][iy]]->Fill(miscalib_[ee],
00335 newCalibs_[ee]);
00336
00337 }
00338 }
00339
00340 EcalCondHeader header;
00341 header.method_="phi symmetry";
00342 header.version_="0";
00343 header.datasource_="testdata";
00344 header.since_=1;
00345 header.tag_="unknown";
00346 header.date_="Mar 24 1973";
00347
00348 EcalIntercalibConstantsXMLTranslator::writeXML(newcalibfile,header,
00349 newCalibs_ );
00350
00351 eehisto.Write();
00352 ebhisto.Write();
00353 ehistof.Close();
00354
00355 fillConstantsHistos();
00356
00357 outResidHistos();
00358
00359 }
00360
00361
00362
00363
00364 void PhiSymmetryCalibration_step2_SM::fillConstantsHistos(){
00365
00366 TFile f("CalibHistos.root","recreate");
00367
00368 TH2F barreletamap("barreletamap","barreletamap",171, -85,86,100,0.,2.);
00369 TH2F barreletamapraw("barreletamapraw","barreletamapraw",171, -85,86,100,0.,2.);
00370
00371 TH2F barrelmapold("barrelmapold","barrelmapold",360,1.,361.,171,-85.,86.);
00372 TH2F barrelmapnew("barrelmapnew","barrelmapnew",360,1.,361.,171,-85.,86.);
00373 TH2F barrelmapratio("barrelmapratio","barrelmapratio",360,1.,361.,171,-85.,86.);
00374
00375 TH1F rawconst_endc_h("rawconst_endc","rawconst_endc",100,0.,2.);
00376 TH1F const_endc_h("const_endc","const_endc",100,0.,2.);
00377
00378 TH1F oldconst_endc_h("oldconst_endc","oldconst_endc;oldCalib;",200,0,2);
00379 TH2F newvsraw_endc_h("newvsraw_endc","newvsraw_endc;rawConst;newCalib",200,0,2,200,0,2);
00380
00381 TH2F endcapmapold_plus("endcapmapold_plus","endcapmapold_plus",100,1.,101.,100,1.,101.);
00382 TH2F endcapmapnew_plus("endcapmapnew_plus","endcapmapnew_plus",100,1.,101.,100,1.,101.);
00383 TH2F endcapmapratio_plus("endcapmapratio_plus","endcapmapratio_plus",100,1.,101.,100,1.,101.);
00384
00385 TH2F endcapmapold_minus("endcapmapold_minus","endcapmapold_minus",100,1.,101.,100,1.,101.);
00386 TH2F endcapmapnew_minus("endcapmapnew_minus","endcapmapnew_minus",100,1.,101.,100,1.,101.);
00387 TH2F endcapmapratio_minus("endcapmapratio_minus","endcapmapratio_minus",100,1.,101.,100,1.,101.);
00388
00389
00390 for (int sign=0; sign<kSides; sign++) {
00391
00392 int thesign = sign==1 ? 1:-1;
00393
00394 for (int ieta=0; ieta<kBarlRings; ieta++) {
00395 for (int iphi=0; iphi<kBarlWedges; iphi++) {
00396 if(e_.goodCell_barl[ieta][iphi][sign]){
00397
00398 EBDetId eb(thesign*( ieta+1 ), iphi+1);
00399
00400
00401 barreletamap.Fill(ieta*thesign + thesign,newCalibs_[eb]);
00402 barreletamapraw.Fill(ieta*thesign + thesign,rawconst_barl[ieta][iphi][sign]);
00403
00404 barrelmapold.Fill(iphi+1,ieta*thesign + thesign, oldCalibs_[eb]);
00405 barrelmapnew.Fill(iphi+1,ieta*thesign + thesign, newCalibs_[eb]);
00406 barrelmapratio.Fill(iphi+1,ieta*thesign + thesign, newCalibs_[eb]/oldCalibs_[eb]);
00407 }
00408 }
00409 }
00410
00411 for (int ix=0; ix<kEndcWedgesX; ix++) {
00412 for (int iy=0; iy<kEndcWedgesY; iy++) {
00413 if (e_.goodCell_endc[ix][iy][sign]){
00414 if (! EEDetId::validDetId(ix+1, iy+1,thesign)) continue;
00415 EEDetId ee(ix+1, iy+1,thesign);
00416
00417 rawconst_endc_h.Fill(rawconst_endc[ix][iy][sign]);
00418 const_endc_h.Fill(newCalibs_[ee]);
00419 oldconst_endc_h.Fill(oldCalibs_[ee]);
00420 newvsraw_endc_h.Fill(rawconst_endc[ix][iy][sign],newCalibs_[ee]);
00421
00422 if(sign==1){
00423 endcapmapold_plus.Fill(ix+1,iy+1,oldCalibs_[ee]);
00424 endcapmapnew_plus.Fill(ix+1,iy+1,newCalibs_[ee]);
00425 endcapmapratio_plus.Fill(ix+1,iy+1,newCalibs_[ee]/oldCalibs_[ee]);
00426 }
00427 else{
00428 endcapmapold_minus.Fill(ix+1,iy+1,oldCalibs_[ee]);
00429 endcapmapnew_minus.Fill(ix+1,iy+1,newCalibs_[ee]);
00430 endcapmapratio_minus.Fill(ix+1,iy+1,newCalibs_[ee]/oldCalibs_[ee]);
00431 }
00432
00433 }
00434 }
00435 }
00436
00437 }
00438
00439 barreletamap.Write();
00440 barreletamapraw.Write();
00441 rawconst_endc_h.Write();
00442 const_endc_h.Write();
00443 oldconst_endc_h.Write();
00444 newvsraw_endc_h.Write();
00445 barrelmapold.Write();
00446 barrelmapnew.Write();
00447 barrelmapratio.Write();
00448 endcapmapold_plus.Write();
00449 endcapmapnew_plus.Write();
00450 endcapmapratio_plus.Write();
00451 endcapmapold_minus.Write();
00452 endcapmapnew_minus.Write();
00453 endcapmapratio_minus.Write();
00454
00455 f.Close();
00456 }
00457
00458
00459
00460
00461
00462
00463
00464
00465
00466
00467
00468 void PhiSymmetryCalibration_step2_SM::fillHistos()
00469 {
00470 TFile f("PhiSymmetryCalibration.root","recreate");
00471
00472 std::vector<TH1F*> etsum_barl_histos(kBarlRings);
00473 std::vector<TH1F*> esum_barl_histos(kBarlRings);
00474
00475
00476 for (int ieta=0; ieta<kBarlRings; ieta++) {
00477 float low=999999.;
00478 float high=0.;
00479 float low_e=999999.;
00480 float high_e=0.;
00481
00482 for (int iphi=0; iphi<kBarlWedges; iphi++) {
00483 for (int sign=0; sign<kSides; sign++) {
00484 float etsum = etsum_barl_[ieta][iphi][sign];
00485 if (etsum<low && etsum!=0.) low=etsum;
00486 if (etsum>high) high=etsum;
00487
00488 float esum = esum_barl_[ieta][iphi][sign];
00489 if (esum<low_e && esum!=0.) low_e=esum;
00490 if (esum>high_e) high_e=esum;
00491 }
00492 }
00493
00494 ostringstream t;
00495 t << "etsum_barl_" << ieta+1;
00496 etsum_barl_histos[ieta]=new TH1F(t.str().c_str(),"",50,low-.2*low,high+.1*high);
00497 t.str("");
00498
00499 t << "esum_barl_" << ieta+1;
00500 esum_barl_histos[ieta]=new TH1F(t.str().c_str(),"",50,low_e-.2*low_e,high_e+.1*high_e);
00501 t.str("");
00502
00503
00504 etsumMean_barl_[ieta]=0.;
00505 esumMean_barl_[ieta]=0.;
00506
00507 for (int iphi=0; iphi<kBarlWedges; iphi++) {
00508 for (int sign=0; sign<kSides; sign++) {
00509
00510
00511 int iphi_r = int(iphi/nscx);
00512
00513 if( !(iphi%nscx)){
00514
00515 etsum_barl_SM_[ieta][iphi_r][sign] = etsum_barl_SM_[ieta][iphi_r][sign]*nscx/(nscx- nBads_barl_SM_[ieta][iphi_r][sign]);
00516
00517
00518 etsumMean_barl_SM_[ieta] +=
00519 etsum_barl_SM_[ieta][iphi_r][sign]/(2*int(kBarlWedges/nscx));
00520 }
00521 if(e_.goodCell_barl[ieta][iphi][sign]){
00522 float etsum = etsum_barl_[ieta][iphi][sign];
00523 float esum = esum_barl_[ieta][iphi][sign];
00524 etsum_barl_histos[ieta]->Fill(etsum);
00525 esum_barl_histos[ieta]->Fill(esum);
00526 etsumMean_barl_[ieta]+=etsum;
00527
00528
00529 esumMean_barl_[ieta]+=esum;
00530 }
00531 }
00532 }
00533
00534 etsum_barl_histos[ieta]->Write();
00535 esum_barl_histos[ieta]->Write();
00536 etsumMean_barl_[ieta]/=(720.-e_.nBads_barl[ieta]);
00537 esumMean_barl_[ieta]/=(720.-e_.nBads_barl[ieta]);
00538
00539 delete etsum_barl_histos[ieta];
00540 delete esum_barl_histos[ieta];
00541 }
00542
00543 std::vector<TH1F*> etsum_endc_histos(kEndcEtaRings);
00544 std::vector<TH1F*> etsum_endc_uncorr_histos(kEndcEtaRings);
00545 std::vector<TH1F*> esum_endc_histos(kEndcEtaRings);
00546
00547 std::vector<TH2F*> etsumvsarea_endc_histos(kEndcEtaRings);
00548 std::vector<TH2F*> esumvsarea_endc_histos(kEndcEtaRings);
00549
00550
00551 for (int ring=0; ring<kEndcEtaRings; ring++) {
00552
00553 float low=FLT_MAX;
00554 float low_uncorr=FLT_MAX;
00555 float high=0.;
00556 float high_uncorr=0;
00557 float low_e=FLT_MAX;
00558 float high_e=0.;
00559 float low_a=1.;
00560 float high_a=0.;
00561 for (int ix=0; ix<kEndcWedgesX; ix++) {
00562 for (int iy=0; iy<kEndcWedgesY; iy++) {
00563 if (e_.endcapRing_[ix][iy]==ring) {
00564 for (int sign=0; sign<kSides; sign++) {
00565 float etsum = etsum_endc_[ix][iy][sign];
00566 if (etsum<low && etsum!=0.) low=etsum;
00567 if (etsum>high) high=etsum;
00568
00569 float etsum_uncorr = etsum_endc_uncorr[ix][iy][sign];
00570 if (etsum_uncorr<low_uncorr && etsum_uncorr!=0.) low_uncorr=etsum_uncorr;
00571 if (etsum_uncorr>high_uncorr) high_uncorr=etsum_uncorr;
00572
00573 float esum = esum_endc_[ix][iy][sign];
00574 if (esum<low_e && esum!=0.) low_e=esum;
00575 if (esum>high_e) high_e=esum;
00576
00577 float area = e_.cellArea_[ix][iy];
00578 if (area<low_a) low_a=area;
00579 if (area>high_a) high_a=area;
00580 }
00581 }
00582 }
00583 }
00584
00585 ostringstream t;
00586 t<<"etsum_endc_" << ring+1;
00587 etsum_endc_histos[ring]= new TH1F(t.str().c_str(),"",50,low-.2*low,high+.1*high);
00588 t.str("");
00589
00590 t<<"etsum_endc_uncorr_" << ring+1;
00591 etsum_endc_uncorr_histos[ring]= new TH1F(t.str().c_str(),"",50,low_uncorr-.2*low_uncorr,high_uncorr+.1*high_uncorr);
00592 t.str("");
00593
00594 t<<"esum_endc_" << ring+1;
00595 esum_endc_histos[ring]= new TH1F(t.str().c_str(),"",50,low_e-.2*low_e,high_e+.1*high_e);
00596 t.str("");
00597
00598 t<<"etsumvsarea_endc_" << ring+1;
00599 etsumvsarea_endc_histos[ring]= new TH2F(t.str().c_str(),";A_{#eta#phi};#Sigma E_{T}",50,low_a,high_a,50,low,high);
00600 t.str("");
00601
00602 t<<"esumvsarea_endc_" << ring+1;
00603 esumvsarea_endc_histos[ring]= new TH2F(t.str().c_str(),";A_{#eta#phi};#Sigma E",50,low_a,high_a,50,low_e,high_e);
00604 t.str("");
00605
00606
00607 etsumMean_endc_[ring]=0.;
00608 esumMean_endc_[ring]=0.;
00609 for (int ix=0; ix<kEndcWedgesX; ix++) {
00610 for (int iy=0; iy<kEndcWedgesY; iy++) {
00611 if (e_.endcapRing_[ix][iy]==ring) {
00612 for (int sign=0; sign<kSides; sign++) {
00613 if(e_.goodCell_endc[ix][iy][sign]){
00614 float etsum = etsum_endc_[ix][iy][sign];
00615 float esum = esum_endc_[ix][iy][sign];
00616 float etsum_uncorr = etsum_endc_uncorr[ix][iy][sign];
00617 etsum_endc_histos[ring]->Fill(etsum);
00618 etsum_endc_uncorr_histos[ring]->Fill(etsum_uncorr);
00619 esum_endc_histos[ring]->Fill(esum);
00620
00621 float area = e_.cellArea_[ix][iy];
00622 etsumvsarea_endc_histos[ring]->Fill(area,etsum);
00623 esumvsarea_endc_histos[ring]->Fill(area,esum);
00624
00625 etsumMean_endc_[ring]+=etsum;
00626 esumMean_endc_[ring]+=esum;
00627 }
00628 }
00629 }
00630 }
00631 }
00632
00633 etsum_endc_histos[ring]->Write();
00634 etsum_endc_uncorr_histos[ring]->Write();
00635 esum_endc_histos[ring]->Write();
00636 etsumMean_endc_[ring]/=(float(e_.nRing_[ring]*2-e_.nBads_endc[ring]));
00637 esumMean_endc_[ring]/=(float(e_.nRing_[ring]*2-e_.nBads_endc[ring]));
00638 etsumvsarea_endc_histos[ring]->Write();
00639 esumvsarea_endc_histos[ring]->Write();
00640
00641 delete etsum_endc_histos[ring];
00642 delete etsum_endc_uncorr_histos[ring];
00643 delete esum_endc_histos[ring];
00644 delete etsumvsarea_endc_histos[ring];
00645 delete esumvsarea_endc_histos[ring];
00646 }
00647
00648
00649
00650 TH2F barreletamap("barreletamap","barreletamap",171, -85,86,100,0,2);
00651 TH2F barrelmap("barrelmap","barrelmap - #frac{#Sigma E_{T}}{<#Sigma E_{T}>_{0}}",360,1,360, 171, -85,86);
00652 TH2F barrelmap_e("barrelmape","barrelmape - #frac{#Sigma E}{<#Sigma E>_{0}}",360,1,360, 171, -85,86);
00653 TH2F barrelmap_divided("barrelmapdiv","barrelmapdivided - #frac{#Sigma E_{T}}{hits}",360,1,360,171,-85,86);
00654 TH2F barrelmap_e_divided("barrelmapediv","barrelmapedivided - #frac{#Sigma E}{hits}",360,1,360,171,-85,86);
00655 TH2F endcmap_plus_corr("endcapmapplus_corrected","endcapmapplus - #frac{#Sigma E_{T}}{<#Sigma E_{T}>_{38}}",100,1,101,100,1,101);
00656 TH2F endcmap_minus_corr("endcapmapminus_corrected","endcapmapminus - #frac{#Sigma E_{T}}{<#Sigma E_{T}>_{38}}",100,1,101,100,1,101);
00657 TH2F endcmap_plus_uncorr("endcapmapplus_uncorrected","endcapmapplus_uncor - #frac{#Sigma E_{T}}{<#Sigma E_{T}>_{38}}",100,1,101,100,1,101);
00658 TH2F endcmap_minus_uncorr("endcapmapminus_uncorrected","endcapmapminus_uncor - #frac{#Sigma E_{T}}{<#Sigma E_{T}>_{38}}",100,1,101,100,1,101);
00659 TH2F endcmap_e_plus("endcapmapeplus","endcapmapeplus - #frac{#Sigma E}{<#Sigma E>_{38}}",100,1,101,100,1,101);
00660 TH2F endcmap_e_minus("endcapmapeminus","endcapmapeminus - #frac{#Sigma E}{<#Sigma E>_{38}}",100,1,101,100,1,101);
00661
00662 for (int sign=0; sign<kSides; sign++) {
00663
00664 int thesign = sign==1 ? 1:-1;
00665
00666 for (int ieta=0; ieta<kBarlRings; ieta++) {
00667 for (int iphi=0; iphi<kBarlWedges; iphi++) {
00668 if(e_.goodCell_barl[ieta][iphi][sign]){
00669 barrelmap.Fill(iphi+1,ieta*thesign + thesign, etsum_barl_[ieta][iphi][sign]/etsumMean_barl_[0]);
00670 barrelmap_e.Fill(iphi+1,ieta*thesign + thesign, esum_barl_[ieta][iphi][sign]/esumMean_barl_[0]);
00671 if (!nhits_barl_[ieta][iphi][sign]) nhits_barl_[ieta][iphi][sign] =1;
00672 barrelmap_divided.Fill( iphi+1,ieta*thesign + thesign, etsum_barl_[ieta][iphi][sign]/nhits_barl_[ieta][iphi][sign]);
00673 barrelmap_e_divided.Fill( iphi+1,ieta*thesign + thesign, esum_barl_[ieta][iphi][sign]/nhits_barl_[ieta][iphi][sign]);
00674
00675
00676 barreletamap.Fill(ieta*thesign + thesign,etsum_barl_[ieta][iphi][sign]/etsumMean_barl_[0]);
00677 }
00678 }
00679 }
00680
00681 for (int ix=0; ix<kEndcWedgesX; ix++) {
00682 for (int iy=0; iy<kEndcWedgesY; iy++) {
00683 if (sign==1) {
00684 endcmap_plus_corr.Fill(ix+1,iy+1,etsum_endc_[ix][iy][sign]/etsumMean_endc_[38]);
00685 endcmap_plus_uncorr.Fill(ix+1,iy+1,etsum_endc_uncorr[ix][iy][sign]/etsumMean_endc_[38]);
00686 endcmap_e_plus.Fill(ix+1,iy+1,esum_endc_[ix][iy][sign]/esumMean_endc_[38]);
00687 }
00688 else{
00689 endcmap_minus_corr.Fill(ix+1,iy+1,etsum_endc_[ix][iy][sign]/etsumMean_endc_[38]);
00690 endcmap_minus_uncorr.Fill(ix+1,iy+1,etsum_endc_uncorr[ix][iy][sign]/etsumMean_endc_[38]);
00691 endcmap_e_minus.Fill(ix+1,iy+1,esum_endc_[ix][iy][sign]/esumMean_endc_[38]);
00692 }
00693 }
00694 }
00695
00696 }
00697
00698
00699
00700 barreletamap.Write();
00701 barrelmap_divided.Write();
00702 barrelmap.Write();
00703 barrelmap_e_divided.Write();
00704 barrelmap_e.Write();
00705 endcmap_plus_corr.Write();
00706 endcmap_minus_corr.Write();
00707 endcmap_plus_uncorr.Write();
00708 endcmap_minus_uncorr.Write();
00709 endcmap_e_plus.Write();
00710 endcmap_e_minus.Write();
00711
00712
00713 vector<TH1F*> etavsphi_endc(kEndcEtaRings);
00714 vector<TH1F*> areavsphi_endc(kEndcEtaRings);
00715 vector<TH1F*> etsumvsphi_endcp_corr(kEndcEtaRings);
00716 vector<TH1F*> etsumvsphi_endcm_corr(kEndcEtaRings);
00717 vector<TH1F*> etsumvsphi_endcp_uncorr(kEndcEtaRings);
00718 vector<TH1F*> etsumvsphi_endcm_uncorr(kEndcEtaRings);
00719 vector<TH1F*> esumvsphi_endcp(kEndcEtaRings);
00720 vector<TH1F*> esumvsphi_endcm(kEndcEtaRings);
00721
00722 std::vector<TH1F*> deltaeta_histos(kEndcEtaRings);
00723 std::vector<TH1F*> deltaphi_histos(kEndcEtaRings);
00724
00725 for(int ring =0; ring<kEndcEtaRings;++ring){
00726
00727 ostringstream t;
00728 t<< "etavsphi_endc_"<<ring;
00729 etavsphi_endc[ring] = new TH1F(t.str().c_str(), t.str().c_str(),e_.nRing_[ring],0,e_.nRing_[ring]);
00730 t.str("");
00731
00732 t<< "areavsphi_endc_"<<ring;
00733 areavsphi_endc[ring] = new TH1F(t.str().c_str(), t.str().c_str(),e_.nRing_[ring],0,e_.nRing_[ring]);
00734 t.str("");
00735
00736 t<< "etsumvsphi_endcp_corr_"<<ring;
00737 etsumvsphi_endcp_corr[ring] = new TH1F(t.str().c_str(), t.str().c_str(),e_.nRing_[ring],0,e_.nRing_[ring]);
00738 t.str("");
00739
00740 t << "etsumvsphi_endcm_corr_"<<ring;
00741 etsumvsphi_endcm_corr[ring] = new TH1F(t.str().c_str(), t.str().c_str(),e_.nRing_[ring],0,e_.nRing_[ring]);
00742 t.str("");
00743
00744 t << "etsumvsphi_endcp_uncorr_"<<ring;
00745 etsumvsphi_endcp_uncorr[ring] = new TH1F(t.str().c_str(), t.str().c_str(),e_.nRing_[ring],0,e_.nRing_[ring]);
00746 t.str("");
00747
00748 t << "etsumvsphi_endcm_uncorr_"<<ring;
00749 etsumvsphi_endcm_uncorr[ring] = new TH1F(t.str().c_str(), t.str().c_str(),e_.nRing_[ring],0,e_.nRing_[ring]);
00750 t.str("");
00751
00752 t << "esumvsphi_endcp_"<<ring;
00753 esumvsphi_endcp[ring] = new TH1F(t.str().c_str(), t.str().c_str(),e_.nRing_[ring],0,e_.nRing_[ring]);
00754 t.str("");
00755
00756 t << "esumvsphi_endcm_"<<ring;
00757 esumvsphi_endcm[ring] = new TH1F(t.str().c_str(), t.str().c_str(),e_.nRing_[ring],0,e_.nRing_[ring]);
00758 t.str("");
00759
00760 t << "deltaeta_" << ring;
00761 deltaeta_histos[ring]= new TH1F(t.str().c_str(),"",50,-.1,.1);
00762 t.str("");
00763 t << "deltaphi_" << ring;
00764 deltaphi_histos[ring]= new TH1F(t.str().c_str(),"",50,-.1,.1);
00765 t.str("");
00766 }
00767
00768 for (int ix=0; ix<kEndcWedgesX; ix++) {
00769 for (int iy=0; iy<kEndcWedgesY; iy++) {
00770
00771 int ring = e_.endcapRing_[ix][iy];
00772 if (ring!=-1) {
00773 int iphi_endc=-1;
00774 for (int ip=0; ip<e_.nRing_[ring]; ip++) {
00775 if (e_.cellPhi_[ix][iy]==e_.phi_endc_[ip][ring]) iphi_endc=ip;
00776 }
00777
00778 if(iphi_endc!=-1){
00779 for (int sign=0; sign<kSides; sign++) {
00780 if(e_.goodCell_endc[ix][iy][sign]){
00781 if (sign==1){
00782 etsumvsphi_endcp_corr[ring]->Fill(iphi_endc,etsum_endc_[ix][iy][sign]);
00783 etsumvsphi_endcp_uncorr[ring]->Fill(iphi_endc,etsum_endc_uncorr[ix][iy][sign]);
00784 esumvsphi_endcp[ring]->Fill(iphi_endc,esum_endc_[ix][iy][sign]);
00785 } else {
00786 etsumvsphi_endcm_corr[ring]->Fill(iphi_endc,etsum_endc_[ix][iy][sign]);
00787 etsumvsphi_endcm_uncorr[ring]->Fill(iphi_endc,etsum_endc_uncorr[ix][iy][sign]);
00788 esumvsphi_endcm[ring]->Fill(iphi_endc,esum_endc_[ix][iy][sign]);
00789 }
00790 }
00791 }
00792 etavsphi_endc[ring]->Fill(iphi_endc,e_.cellPos_[ix][iy].eta());
00793 areavsphi_endc[ring]->Fill(iphi_endc,e_.cellArea_[ix][iy]);
00794 }
00795
00796 }
00797 }
00798 }
00799
00800
00801
00802 for(int ring =0; ring<kEndcEtaRings;++ring){
00803
00804 etavsphi_endc[ring]->Write();
00805 areavsphi_endc[ring]->Write();
00806 etsumvsphi_endcp_corr[ring]->Write();
00807 etsumvsphi_endcm_corr[ring]->Write();
00808 etsumvsphi_endcp_uncorr[ring]->Write();
00809 etsumvsphi_endcm_uncorr[ring]->Write();
00810 esumvsphi_endcp[ring]->Write();
00811 esumvsphi_endcm[ring]->Write();
00812 deltaeta_histos[ring]->Write();
00813 deltaphi_histos[ring]->Write();
00814
00815
00816 delete etsumvsphi_endcp_corr[ring];
00817 delete etsumvsphi_endcm_corr[ring];
00818 delete etsumvsphi_endcp_uncorr[ring];
00819 delete etsumvsphi_endcm_uncorr[ring];
00820 delete etavsphi_endc[ring];
00821 delete areavsphi_endc[ring];
00822 delete esumvsphi_endcp[ring];
00823 delete esumvsphi_endcm[ring];
00824 delete deltaeta_histos[ring];
00825 delete deltaphi_histos[ring];
00826 }
00827
00828
00829 f.Close();
00830 }
00831
00832
00833 void PhiSymmetryCalibration_step2_SM::readEtSums(){
00834
00835
00836
00837
00838 int ieta,iphi,sign,ix,iy,dummy;
00839 double etsum;
00840 unsigned int nhits;
00841 std::ifstream etsum_barl_in("etsum_barl.dat", ios::in);
00842 while ( etsum_barl_in >> dummy >> ieta >> iphi >> sign >> etsum >> nhits ) {
00843 etsum_barl_[ieta][iphi][sign]+=etsum;
00844 nhits_barl_[ieta][iphi][sign]+=nhits;
00845
00846
00847 int iphi_r = int(iphi/nscx);
00848 etsum_barl_SM_[ieta][iphi_r][sign]+= etsum;
00849
00850
00851
00852
00853 }
00854
00855 std::ifstream etsum_endc_in("etsum_endc.dat", ios::in);
00856 while ( etsum_endc_in >> dummy >> ix >> iy >> sign >> etsum >> nhits>>dummy ) {
00857 etsum_endc_[ix][iy][sign]+=etsum;
00858 nhits_endc_[ix][iy][sign]+=nhits;
00859
00860 }
00861
00862 std::ifstream k_barl_in("k_barl.dat", ios::in);
00863 for (int ieta=0; ieta<kBarlRings; ieta++) {
00864 k_barl_in >> dummy >> k_barl_[ieta];
00865 }
00866
00867 std::ifstream k_endc_in("k_endc.dat", ios::in);
00868 for (int ring=0; ring<kEndcEtaRings; ring++) {
00869 k_endc_in >> dummy >> k_endc_[ring];
00870 }
00871
00872
00873 }
00874
00875
00876
00877 void PhiSymmetryCalibration_step2_SM::setupResidHistos(){
00878
00879 miscal_resid_barl_histos.resize(kBarlRings);
00880 correl_barl_histos.resize(kBarlRings);
00881
00882
00883 for (int ieta=0; ieta<kBarlRings; ieta++) {
00884 ostringstream t1;
00885 t1<<"mr_barl_"<<ieta+1;
00886 miscal_resid_barl_histos[ieta] = new TH1F(t1.str().c_str(),"",100,0.,2.);
00887 ostringstream t2;
00888 t2<<"co_barl_"<<ieta+1;
00889 correl_barl_histos[ieta] = new TH2F(t2.str().c_str(),"",50,.5,1.5,50,.5,1.5);
00890 }
00891
00892 miscal_resid_endc_histos.resize(kEndcEtaRings);
00893 correl_endc_histos.resize(kEndcEtaRings);
00894
00895
00896 for (int ring=0; ring<kEndcEtaRings; ring++) {
00897 ostringstream t1;
00898 t1<<"mr_endc_"<< ring+1;
00899 miscal_resid_endc_histos[ring] = new TH1F(t1.str().c_str(),"",100,0.,2.);
00900 ostringstream t2;
00901 t2<<"co_endc_"<<ring+1;
00902 correl_endc_histos[ring] = new TH2F(t2.str().c_str(),"",50,.5,1.5,50,.5,1.5);
00903 }
00904
00905
00906
00907
00908
00909 }
00910
00911
00912 void PhiSymmetryCalibration_step2_SM::outResidHistos(){
00913
00914
00915 TFile f("PhiSymmetryCalibration_miscal_resid.root","recreate");
00916 for (int ieta=0; ieta<85; ieta++) {
00917 miscal_resid_barl_histos[ieta]->Write();
00918 correl_barl_histos[ieta]->Write();
00919
00920 delete miscal_resid_barl_histos[ieta];
00921 delete correl_barl_histos[ieta];
00922 }
00923
00924 for (int ring=0; ring<39; ring++) {
00925 miscal_resid_endc_histos[ring]->Write();
00926 correl_endc_histos[ring]->Write();
00927
00928 delete miscal_resid_endc_histos[ring];
00929 delete correl_endc_histos[ring];
00930
00931 }
00932 f.Close();
00933 }