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