CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/Calibration/HcalCalibAlgos/src/hcalCalib.cc

Go to the documentation of this file.
00001 #define hcalCalib_cxx
00002 
00003 //  TSelector-based code for getting the HCAL resp. correction
00004 //  from physics events. Works for DiJet and IsoTrack calibration.
00005 //
00006 //  Anton Anastassov (Northwestern)
00007 //  Email: aa@fnal.gov
00008 //
00009 // $Id: hcalCalib.cc,v 1.8 2010/10/03 20:38:58 elmer Exp $
00010 //
00011 
00012 #include "Calibration/HcalCalibAlgos/interface/hcalCalib.h"
00013 #include "Calibration/HcalCalibAlgos/interface/hcalCalibUtils.h"
00014 
00015 #include <TH2.h>
00016 #include <TStyle.h>
00017 #include "TFile.h"
00018 
00019 #include <iostream>
00020 #include <fstream>
00021 
00022 #include <sstream>
00023 #include <string>
00024 
00025 #include <map>
00026 #include <numeric>
00027 #include <algorithm>
00028 #include <set>
00029 
00030 #include "Calibration/Tools/interface/MinL3AlgoUniv.h"
00031 
00032 #include "TMatrixF.h"
00033 #include "TMatrixD.h"
00034 #include "TDecompSVD.h"
00035 #include "TDecompQRH.h"
00036 
00037 #include "DataFormats/DetId/interface/DetId.h"
00038 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00039 
00040 
00041 using namespace std;
00042 
00043 UInt_t nEvents;
00044 
00045 TFile* histoFile;
00046 
00047 // sanity check histograms
00048 TH1F* h1_trkP;
00049 TH1F* h1_allTrkP;
00050 
00051 TH1F* h1_selTrkP_iEta10;
00052 
00053 TH1F* h1_rawSumE;
00054 TH1F* h1_rawResp;
00055 TH1F* h1_corResp;
00056 TH1F* h1_rawRespBarrel;
00057 TH1F* h1_corRespBarrel;
00058 TH1F* h1_rawRespEndcap;
00059 TH1F* h1_corRespEndcap;
00060 TH1F* h1_numEventsTwrIEta;
00061 
00062 TH2F* h2_dHitRefBarrel;
00063 TH2F* h2_dHitRefEndcap;
00064 
00065 // histograms based on iEta, iPhi of refPosition forthe cluster (at the moment: hottest tower)
00066 // expect range |iEta|<=24 (to do: add flexibility for arbitrary range)
00067 TH1F* h1_corRespIEta[48];
00068 
00069 
00070 void hcalCalib::Begin(TTree * /*tree*/) {
00071   TString option = GetOption();
00072 
00073   nEvents = 0;
00074 
00075   if (APPLY_PHI_SYM_COR_FLAG && !ReadPhiSymCor()) {
00076     cout << "\nERROR: Failed to read the phi symmetry corrections." << endl;
00077     cout << "Check if the filename is correct. If the corrections are not needed, set the corresponding flag to \"false\"\n" << endl;
00078 
00079     cout << "\nThe program will be terminated\n" << endl;
00080 
00081     exit(1);
00082 
00083   }
00084 
00085     
00086   //  cellEnergies.reserve(1000000);
00087   //  cellIds.reserve(1000000);
00088   //  targetEnergies.reserve(1000000);
00089     
00090       histoFile = new TFile(HISTO_FILENAME.Data(), "RECREATE");
00091 
00092 
00093   h1_trkP    = new TH1F("h1_trkP", "Track momenta; p_{trk} (GeV); Number of tracks", 100, 0, 200);
00094   h1_allTrkP = new TH1F("h1_allTrkP", "Track momenta - all tracks; p_{trk} (GeV); Number of tracks", 100, 0, 200);
00095 
00096   h1_selTrkP_iEta10 = new TH1F("h1_selTrkP_iEta10", "Track momenta - tracks with |iEta|<10; p_{trk} (GeV); Number of tracks", 100, 0, 200); 
00097 
00098 
00099 
00100   if (CALIB_TYPE=="ISO_TRACK") 
00101     h1_rawSumE = new TH1F("h1_rawSumE", "Cluster Energy; E_{cl} (GeV); Number of tracks", 100, 0, 200);
00102   else 
00103     h1_rawSumE = new TH1F("h1_rawSumE", "Cluster Energy; E_{cl} (GeV); Number of tracks", 1000, 0, 2000);   
00104 
00105 
00106 
00107   h1_rawResp = new TH1F("h1_rawResp", "Uncorrected response: |iEta|<24;  E_{had}/p; Number of tracks", 300, 0, 3);
00108   h1_corResp = new TH1F("h1_corResp", "Corrected response: |iEta|<24; E_{had}/p; Number of tracks", 300, 0, 3);
00109     
00110   h1_rawRespBarrel = new TH1F("h1_rawRespBarrel", "Uncorrected response: |iEta|<15;  E_{had}/p; Number of tracks", 300, 0, 3);
00111   h1_corRespBarrel = new TH1F("h1_corRespBarrel", "Corrected response: |iEta|<15; E_{had}/p; Number of tracks", 300, 0, 3);
00112     
00113   h1_rawRespEndcap = new TH1F("h1_rawRespEndcap", "Uncorrected response:  17<|iEta|<24;  E_{had}/p; Number of tracks", 300, 0, 3);
00114   h1_corRespEndcap = new TH1F("h1_corRespEndcap", "Corrected response: 17<|iEta|<24; E_{had}/p; Number of tracks", 300, 0, 3);
00115   
00116   h1_numEventsTwrIEta = new TH1F("h1_numEventsTwrIEta", "h1_numEventsTwrIEta", 80, -40, 40);
00117 
00118   h2_dHitRefBarrel     = new TH2F("h2_dHitRefBarrel","{#Delta}i{#phi} vs {#Delta}i{#eta} of hit and most energetic tower(|i{#eta}|<16);{#Delta}i{#eta}; {#Delta}i{#phi}", 10, -5, 5, 10, -5, 5);
00119   h2_dHitRefEndcap     = new TH2F("h2_dHitRefEndcap","{#Delta}i{#phi} vs {#Delta}i{#eta} of hit and most energetic tower (16<|i{#eta}|<25) ;{#Delta}i{#eta}; {#Delta}i{#phi}", 10, -5, 5, 10, -5, 5);
00120       
00121 
00122 
00123   TString histoName = "isoTrack_";
00124  
00125   for (Int_t i=0; i<48; ++i) {
00126     Long_t iEta; 
00127     if   (i<24) iEta = i-24;
00128     else iEta = i-23;
00129     TString hn = histoName + iEta;
00130     h1_corRespIEta[i] = new TH1F(hn, hn, 300, 0, 3.0);
00131   }
00132 
00133 
00134 }  // end of Begin()
00135 
00136 
00137 
00138 //void hcalCalib::SlaveBegin(TTree * /*tree*/) {
00139 //  TString option = GetOption();
00140 //}
00141 
00142 
00143 Bool_t hcalCalib::Process(Long64_t entry) {
00144  
00145   //  fChain->GetTree()->GetEntry(entry);
00146   GetEntry(entry);
00147 
00148 
00149   set<UInt_t>  uniqueIds; // for testing: check if there are duplicate cells   (AA)
00150 
00151  
00152   Bool_t acceptEvent = kTRUE;
00153    
00154   ++nEvents;
00155 
00156   if (!(nEvents%100000)) cout << "event: " << nEvents << endl;
00157 
00158 
00159   h1_allTrkP->Fill(targetE);
00160    
00161   if (targetE < MIN_TARGET_E || targetE > MAX_TARGET_E) return kFALSE;;
00162 
00163  
00164   // make local copy as the cells may be modified  due to phi/depth sum, phi corrections etc  
00165   vector<TCell> selectCells;
00166 
00167 
00168   if (cells->GetSize()==0) return kFALSE;
00169 
00170   if (CALIB_TYPE=="DI_JET" && probeJetEmFrac > 0.999) return kTRUE;
00171   
00172 
00173   for (Int_t i=0; i<cells->GetSize(); ++i) {
00174     TCell* thisCell =  (TCell*) cells->At(i);
00175 
00176     if (HcalDetId(thisCell->id()).subdet()== HcalOuter)  continue;           // reject HO, make a switch!
00177 
00178     if (HcalDetId(thisCell->id()).subdet() != HcalBarrel &&
00179         HcalDetId(thisCell->id()).subdet() != HcalEndcap &&
00180         HcalDetId(thisCell->id()).subdet() != HcalForward) {
00181       
00182       cout << "Unknown or wrong hcal subdetector: " << HcalDetId(thisCell->id()).subdet() << endl;
00183 
00184     }
00185 
00186 
00187     // Apply phi symmetry corrections if the flag is set
00188     if (APPLY_PHI_SYM_COR_FLAG) thisCell->SetE(phiSymCor[thisCell->id()] * thisCell->e());
00189 
00190     if (thisCell->e() > MIN_CELL_E) selectCells.push_back(*thisCell);
00191   }
00192 
00193 
00194   if (selectCells.size()==0) {
00195     cout << "NO CELLS ABOVE THRESHOLD FOUND FOR TARGET!!!" << endl;
00196   }
00197 
00198 
00199   if (SUM_DEPTHS)  sumDepths(selectCells);
00200   else if (SUM_SMALL_DEPTHS) sumSmallDepths(selectCells); // depth 1,2 in twrs 15,16
00201 
00202 
00203 
00204   // most energetic tower (IsoTracks) or centroid of probe jet (DiJets)
00205   pair<Int_t, UInt_t> refPos;
00206 
00207   Int_t dEtaHitRef = 999;
00208   Int_t dPhiHitRef = 999;
00209 
00210   if (CALIB_TYPE=="ISO_TRACK") { 
00211 
00212     Int_t  iEtaMaxE;  // filled by reference in getIEtaIPhiForHighestE
00213     UInt_t iPhiMaxE;  // 
00214       
00215     getIEtaIPhiForHighestE(selectCells, iEtaMaxE, iPhiMaxE);   
00216  
00217     dEtaHitRef = iEtaMaxE - iEtaHit;
00218     dPhiHitRef = iPhiMaxE - iPhiHit;
00219 
00220     if (dPhiHitRef < -36) dPhiHitRef += 72;
00221     if (dPhiHitRef >  36) dPhiHitRef -= 72;
00222 
00223     if  (iEtaHit*iEtaMaxE < 0) {
00224       if (dEtaHitRef<0) dEtaHitRef += 1;
00225       if (dEtaHitRef>0) dEtaHitRef -= 1;     
00226     }
00227 
00228       
00229     if (abs(iEtaHit)<16) h2_dHitRefBarrel->Fill(dEtaHitRef, dPhiHitRef);
00230     if (abs(iEtaHit)>16 && abs(iEtaHit)<25) h2_dHitRefEndcap->Fill(dEtaHitRef, dPhiHitRef);
00231 
00232     // --------------------------------------------------
00233     // Choice of cluster definition 
00234     //
00235     // fixed size NxN clusters as specified in to config file 
00236     if (!USE_CONE_CLUSTERING) {
00237       if (abs(iEtaMaxE)<16 && HB_CLUSTER_SIZE==3) filterCells3x3(selectCells, iEtaMaxE, iPhiMaxE);
00238       if (abs(iEtaMaxE)>15 && HE_CLUSTER_SIZE==3) filterCells3x3(selectCells, iEtaMaxE, iPhiMaxE);
00239       
00240       if (abs(iEtaMaxE)<16 && HB_CLUSTER_SIZE==5) filterCells5x5(selectCells, iEtaMaxE, iPhiMaxE);
00241       if (abs(iEtaMaxE)>15 && HE_CLUSTER_SIZE==5) filterCells5x5(selectCells, iEtaMaxE, iPhiMaxE);    
00242     }      
00243     else {
00244       //  calculate distance at hcal surface
00245       const GlobalPoint hitPositionHcal(xTrkHcal, yTrkHcal, zTrkHcal); 
00246       filterCellsInCone(selectCells, hitPositionHcal, MAX_CONE_DIST, theCaloGeometry);
00247     }
00248 
00249 
00250 
00251     refPos.first  = iEtaMaxE;
00252     refPos.second = iPhiMaxE;
00253 
00254   }
00255   else if (CALIB_TYPE=="DI_JET") {          // Apply selection cuts on DiJet events here
00256 
00257     if (etVetoJet>MAX_ET_THIRD_JET) acceptEvent = kFALSE;
00258     
00259     Float_t jetsDPhi = probeJetP4->DeltaPhi(*tagJetP4); 
00260     if (fabs(jetsDPhi * 180.0/M_PI) < MIN_DPHI_DIJETS) acceptEvent = kFALSE;
00261 
00262     if (probeJetEmFrac>MAX_PROBEJET_EMFRAC)          acceptEvent = kFALSE;
00263     if (fabs(probeJetP4->Eta())<MIN_PROBEJET_ABSETA) acceptEvent = kFALSE;
00264     if (fabs(tagJetP4->Eta())>MAX_TAGJET_ABSETA)     acceptEvent = kFALSE;
00265     if (fabs(tagJetP4->Et())< MIN_TAGJET_ET)         acceptEvent = kFALSE;
00266 
00267     if (acceptEvent) {
00268 
00269       Int_t  iEtaMaxE;  // filled by reference in getIEtaIPhiForHighestE
00270       UInt_t iPhiMaxE;  // 
00271       
00272       getIEtaIPhiForHighestE(selectCells, iEtaMaxE, iPhiMaxE);   
00273       
00274       // The ref position for the jet is not used in the minimization at this time.
00275       // It will be needed if we attempt to do matrix inversion: then the question is
00276       // which value is better suited: the centroid of the jet or the hottest tower...
00277 
00278       //    refPos.first  = iEtaHit;
00279       //    refPos.second = iPhiHit;
00280       
00281       refPos.first  = iEtaMaxE;
00282       refPos.second = iPhiMaxE;  
00283       
00284       if (abs(iEtaMaxE)>40) acceptEvent = kFALSE;             // for testing :  set as parameter (AA)
00285 
00286     }
00287 
00288   }
00289 
00290 
00291   if (COMBINE_PHI) combinePhi(selectCells);
00292     
00293   // fill the containers for the minimization prcedures
00294   vector<Float_t> energies;
00295   vector<UInt_t>  ids;
00296     
00297   for (vector<TCell>::iterator i_it = selectCells.begin(); i_it!=selectCells.end(); ++i_it) {
00298 
00299     // for testing : fill only unique id's
00300 
00301     if (uniqueIds.insert(i_it->id()).second) {
00302       energies.push_back(i_it->e());
00303       ids.push_back(i_it->id());
00304     }
00305 
00306   }
00307 
00308   if (CALIB_TYPE=="ISO_TRACK") {
00309     if (accumulate(energies.begin(), energies.end(), 0.0) / targetE < MIN_EOVERP) acceptEvent=kFALSE;
00310     if (accumulate(energies.begin(), energies.end(), 0.0) / targetE > MAX_EOVERP) acceptEvent=kFALSE;
00311 
00312     if (emEnergy > MAX_TRK_EME) acceptEvent=kFALSE;
00313 
00314     if (abs(dEtaHitRef)>1 || abs(dPhiHitRef)>1) acceptEvent=kFALSE;
00315 
00316     // Have to check if for |iEta|>20 (and neighboring region) the dPhiHitRef
00317     // should be relaxed to 2. The neighboring towers have dPhi=2...
00318 
00319 
00320   }
00321 
00322     
00323   h1_rawSumE->Fill(accumulate(energies.begin(), energies.end(), 0.0));
00324     
00325      
00326   // here we fill the information for the minimization procedure
00327   if (acceptEvent) { 
00328     cellEnergies.push_back(energies);
00329     cellIds.push_back(ids);
00330     targetEnergies.push_back(targetE);
00331     refIEtaIPhi.push_back(refPos);
00332 
00333     if (abs(refPos.first)<=10) 
00334       h1_selTrkP_iEta10->Fill(targetE);
00335 
00336 
00337   }
00338 
00339 
00340 
00341     
00342   // Clean up
00343   energies.clear();
00344   ids.clear();
00345   selectCells.clear();
00346     
00347   return kTRUE;
00348 }
00349 
00350 //void hcalCalib::SlaveTerminate() {}
00351 
00352 void hcalCalib::Terminate() {
00353 
00354   cout << "\n\nFinished reading the events.\n";
00355   cout << "Number of input objects: " << cellIds.size()  << endl; 
00356   cout << "Performing minimization: depending on selected method can take some time...\n\n";
00357   
00358   for (vector<pair<Int_t, UInt_t> >::iterator it_rp=refIEtaIPhi.begin(); it_rp!=refIEtaIPhi.end(); ++it_rp ) {
00359     Float_t weight = (abs(it_rp->first)<21)? 1.0/72.0 : 1.0/36.0;
00360     h1_numEventsTwrIEta->Fill(it_rp->first, weight);
00361   }
00362  
00363   if (CALIB_METHOD=="MATRIX_INV_OF_ETA_AVE") {
00364     GetCoefFromMtrxInvOfAve();
00365    }    
00366   else if (CALIB_METHOD=="L3" || CALIB_METHOD=="L3_AND_MTRX_INV") {
00367     int eventWeight = 2;  // 2 is the default (try at some point 0,1,2,3)
00368     MinL3AlgoUniv<UInt_t>* thisL3Algo = new MinL3AlgoUniv<UInt_t>(eventWeight);
00369     int numIterations = 10;  // default 10
00370         
00371     solution = thisL3Algo->iterate(cellEnergies, cellIds, targetEnergies, numIterations);
00372 
00373     // in order to handle the case where sumDepths="false", but the flag to sum depths 1,2 in HB towers 15, 16
00374     // is set (sumSmallDepths) we create entries in "solution" to create equal correction here
00375     // for each encountered coef in depth one.
00376 
00377     if (!SUM_DEPTHS && SUM_SMALL_DEPTHS) {
00378       vector<UInt_t> idForSummedCells;
00379 
00380       for (map<UInt_t, Float_t>::iterator m_it = solution.begin(); m_it != solution.end(); ++m_it) {
00381         if (HcalDetId(m_it->first).ietaAbs()!=15 && HcalDetId(m_it->first).ietaAbs()!=16) continue;
00382         if (HcalDetId(m_it->first).subdet()!=HcalBarrel) continue;
00383         if (HcalDetId(m_it->first).depth()==1)
00384           idForSummedCells.push_back(HcalDetId(m_it->first));
00385       }
00386 
00387       for (vector<UInt_t>::iterator v_it=idForSummedCells.begin(); v_it!=idForSummedCells.end(); ++v_it) {
00388         UInt_t addCoefId = HcalDetId(HcalBarrel, HcalDetId(*v_it).ieta(), HcalDetId(*v_it).iphi(), 2);
00389         solution[addCoefId] = solution[*v_it];
00390       }
00391      
00392     }  // end of special treatment for "sumSmallDepths" mode
00393 
00394 
00395     if (CALIB_METHOD=="L3_AND_MTRX_INV") {
00396       GetCoefFromMtrxInvOfAve();      
00397       
00398       // loop over the solution from L3 and multiply by the additional factor from
00399       // the matrix inversion. Set coef outside of the valid calibration region =1.       
00400       for (map<UInt_t, Float_t>::iterator it_s=solution.begin(); it_s!=solution.end(); ++it_s) {
00401         Int_t iEtaSol = HcalDetId(it_s->first).ieta();
00402         if (abs(iEtaSol)<CALIB_ABS_IETA_MIN || abs(iEtaSol)> CALIB_ABS_IETA_MAX) it_s->second = 1.0;
00403         else  it_s->second *= iEtaCoefMap[iEtaSol]; 
00404       }    
00405       
00406     } // if CALIB_METHOD=="L3_AND_MTRX_INV"
00407     
00408   }  // end of L3 or L3 + mtrx inversion minimization
00409 
00410   // done getting the constants -> write the formatted file
00411   makeTextFile();
00412   
00413   // fill some histograms
00414   Float_t rawResp = 0;
00415   Float_t corResp = 0;
00416   Int_t maxIEta = 0;
00417   Int_t minIEta = 999;
00418 
00419 
00420   for (unsigned int i=0; i<cellEnergies.size(); ++i) {
00421     Int_t iEta;    
00422     for (unsigned int j=0; j<(cellEnergies[i]).size(); ++j) {
00423       iEta = HcalDetId(cellIds[i][j]).ieta();
00424       rawResp += (cellEnergies[i])[j];
00425       
00426       if (CALIB_METHOD=="L3_AND_MTRX_INV") {
00427         corResp += (solution[cellIds[i][j]] * cellEnergies[i][j]);
00428       }
00429       else if (CALIB_METHOD=="L3") {
00430         corResp += (solution[cellIds[i][j]] * (cellEnergies[i][j]));
00431       }
00432       else if (CALIB_METHOD=="MATRIX_INV_OF_ETA_AVE") {
00433         corResp += (iEtaCoefMap[iEta]* cellEnergies[i][j]);     
00434       }
00435       
00436       if (maxIEta < abs(iEta) ) maxIEta = abs(iEta);
00437       if (minIEta > abs(iEta) ) minIEta = abs(iEta); 
00438     }
00439 
00440     rawResp /= targetEnergies[i];
00441     corResp /= targetEnergies[i];
00442 
00443 
00444     // fill histograms based on iEta on ref point of the cluster (for now: hottest tower) 
00445     // expect range |iEta|<=24 (to do: add flexibility for arbitrary range)
00446  
00447     if (CALIB_TYPE=="ISO_TRACK") {
00448       Int_t ind = refIEtaIPhi[i].first;  
00449       (ind<0) ? (ind +=24) : (ind +=23);
00450       if (ind>=0 && ind<48) {
00451         h1_corRespIEta[ind]->Fill(corResp);  
00452       }
00453      
00454       // fill histograms for cases where all towers are in the barrel or endcap
00455       if (maxIEta<25) {
00456         h1_rawResp->Fill(rawResp);
00457         h1_corResp->Fill(corResp);
00458       }
00459       if (maxIEta<15) {
00460         h1_rawRespBarrel->Fill(rawResp);
00461         h1_corRespBarrel->Fill(corResp);
00462       }
00463       else if (maxIEta<25 && minIEta>16) {
00464         h1_rawRespEndcap->Fill(rawResp);
00465         h1_corRespEndcap->Fill(corResp);
00466       }
00467     }  // histograms for isotrack calibration
00468 
00469 
00470     else {
00471      
00472       // put jet plots here
00473 
00474     }
00475 
00476 
00477 
00478     rawResp = 0;
00479     corResp = 0;
00480     maxIEta = 0;
00481     minIEta = 999;
00482   }
00483    
00484 
00485   // save the histograms 
00486   h1_trkP->Write();
00487   h1_allTrkP->Write();
00488 
00489 h1_selTrkP_iEta10->Write();
00490 
00491   h1_rawSumE->Write();
00492   h1_rawResp->Write();
00493   h1_corResp->Write();
00494   h1_rawRespBarrel->Write();
00495   h1_corRespBarrel->Write();
00496   h1_rawRespEndcap->Write();
00497   h1_corRespEndcap->Write();
00498   h1_numEventsTwrIEta->Write();
00499   h2_dHitRefBarrel->Write();
00500   h2_dHitRefEndcap->Write();
00501   for (Int_t i=0; i<48; ++i) {
00502     h1_corRespIEta[i]->Write();
00503   }
00504 
00505 
00506 
00507 
00508 
00509   histoFile->Write();
00510   histoFile->Close();
00511 
00512   cout << "\n Finished calibration.\n " << endl; 
00513 
00514     
00515 }  // end of Terminate()
00516 
00517 
00518 
00519 
00520 //**************************************************************************************************************
00521 
00522 void hcalCalib::GetCoefFromMtrxInvOfAve() {
00523 
00524     // these maps are keyed by iEta        
00525     map<Int_t, Float_t> aveTargetE;     
00526     map<Int_t, Int_t>   nEntries;     // count hits
00527         
00528     //  iEtaRef  iEtaCell, energy 
00529     map<Int_t, map<Int_t, Float_t> > aveHitE; // add energies in the loop, normalize after that
00530  
00531     for (unsigned int i=0; i<cellEnergies.size(); ++i) {
00532       Int_t iEtaRef = refIEtaIPhi[i].first;
00533       aveTargetE[iEtaRef] += targetEnergies[i];
00534       nEntries[iEtaRef]++;
00535 
00536       // for hybrid method: matrix inv of averages preceeded by L3
00537       if (CALIB_METHOD=="L3_AND_MTRX_INV") {
00538         for (unsigned int j=0; j<(cellEnergies[i]).size(); ++j) {
00539           aveHitE[iEtaRef][HcalDetId(cellIds[i][j]).ieta()] += (solution[cellIds[i][j]] * cellEnergies[i][j]);
00540         }          
00541       }
00542       else if (CALIB_METHOD=="MATRIX_INV_OF_ETA_AVE") {
00543         for (unsigned int j=0; j<(cellEnergies[i]).size(); ++j) {
00544           aveHitE[iEtaRef][HcalDetId(cellIds[i][j]).ieta()] += cellEnergies[i][j];
00545         }          
00546       }
00547 
00548     } // end of loop of entries
00549 
00550            
00551     // scale by number of entries to get the averages
00552     Float_t norm = 1.0;
00553     for (map<Int_t, Float_t>::iterator m_it = aveTargetE.begin(); m_it!=aveTargetE.end(); ++m_it){            
00554       Int_t iEta = m_it->first;
00555       norm = (nEntries[iEta] > 0)? 1.0/(nEntries[iEta]) : 1.0;
00556       aveTargetE[iEta] *= norm;
00557 
00558       map<Int_t, Float_t>::iterator n_it = (aveHitE[iEta]).begin();
00559 
00560       Float_t sumRawE = 0;
00561       for (; n_it!=(aveHitE[iEta]).end(); ++n_it) {
00562         (n_it->second) *= norm;
00563         sumRawE += (n_it->second);
00564       }
00565         
00566     }  // end of scaling by number of entries
00567   
00568     Int_t ONE_SIDE_IETA_RANGE = CALIB_ABS_IETA_MAX - CALIB_ABS_IETA_MIN +1;
00569 
00570     // conversion from iEta to index for the linear system
00571     // contains elements only in the valid range for *matrix inversion*
00572     vector<Int_t> iEtaList;
00573 
00574     for (Int_t i=-CALIB_ABS_IETA_MAX; i<=CALIB_ABS_IETA_MAX; ++i) {
00575       if (abs(i)<CALIB_ABS_IETA_MIN) continue;
00576       iEtaList.push_back(i); 
00577     }   
00578 
00579     TMatrixD A(2*ONE_SIDE_IETA_RANGE, 2*ONE_SIDE_IETA_RANGE);
00580     TMatrixD b(2*ONE_SIDE_IETA_RANGE, 1);
00581     TMatrixD x(2*ONE_SIDE_IETA_RANGE, 1);
00582         
00583     for (Int_t i=0; i<2*ONE_SIDE_IETA_RANGE; ++i) {
00584       for (Int_t j=0; j<2*ONE_SIDE_IETA_RANGE; ++j) {
00585         A(i, j) = 0.0;
00586       }
00587     }
00588 
00589     for (UInt_t i=0; i<iEtaList.size(); ++i) {     
00590       b(i, 0) = aveTargetE[iEtaList[i]];        
00591 
00592       map<Int_t, Float_t>::iterator n_it = aveHitE[iEtaList[i]].begin();
00593       for (; n_it!=aveHitE[iEtaList[i]].end(); ++n_it){
00594 
00595         if (fabs(n_it->first)>CALIB_ABS_IETA_MAX || fabs(n_it->first)<CALIB_ABS_IETA_MIN) continue;
00596         Int_t j =  Int_t(find(iEtaList.begin(),iEtaList.end(), n_it->first) - iEtaList.begin());        
00597         A(i, j) = n_it->second;
00598 
00599       }
00600 
00601     }   
00602 
00603     TMatrixD coef = b;
00604     TDecompQRH qrh(A);
00605     Bool_t hasSolution = qrh.MultiSolve(coef);
00606    
00607     if (hasSolution && (CALIB_METHOD=="L3_AND_MTRX_INV" || CALIB_METHOD=="MATRIX_INV_OF_ETA_AVE")) {
00608       for (UInt_t i=0; i<iEtaList.size(); ++i) {
00609         iEtaCoefMap[iEtaList[i]] = coef(i, 0);
00610       }
00611     }
00612 }
00613 
00614 
00615 Bool_t hcalCalib::ReadPhiSymCor() {
00616 
00617   ifstream phiSymFile(PHI_SYM_COR_FILENAME.Data());
00618 
00619   if (!phiSymFile) {
00620     cout << "\nERROR: Can not find file with phi symmetry constants \"" <<  PHI_SYM_COR_FILENAME.Data() << "\"" << endl;
00621     return kFALSE;
00622   }
00623 
00624   // assumes the format used in CSA08, first line is a comment
00625 
00626   Int_t   iEta;
00627   UInt_t  iPhi;
00628   UInt_t  depth;
00629   TString sdName;
00630   UInt_t  detId;
00631 
00632   Float_t value;
00633   HcalSubdetector sd;
00634 
00635   std::string line;
00636 
00637   while (getline(phiSymFile, line)) {
00638 
00639     if(!line.size() || line[0]=='#') continue;
00640 
00641     std::istringstream linestream(line);
00642     linestream >> iEta >> iPhi >> depth >> sdName >> value >> hex >> detId;
00643 
00644     if (sdName=="HB") sd = HcalBarrel;
00645     else if (sdName=="HE") sd = HcalEndcap;
00646     else if (sdName=="HO") sd = HcalOuter;
00647     else if (sdName=="HF") sd = HcalForward;
00648     else {
00649       cout << "\nInvalid detector name in phi symmetry constants file: " << sdName.Data() << endl;
00650       cout << "Check file and rerun!\n" << endl;
00651       return kFALSE;
00652     }
00653 
00654     // check if the data is consistent    
00655   
00656     if (HcalDetId(sd, iEta, iPhi, depth) != HcalDetId(detId)) {
00657       cout << "\nInconsistent info in phi symmetry file: subdet, iEta, iPhi, depth do not match rawId!\n" << endl;
00658       return kFALSE;    
00659     }
00660     if (!HcalDetId::validDetId(sd, iEta, iPhi, depth)) {
00661       cout << "\nInvalid DetId from: iEta=" << iEta << " iPhi=" << iPhi << " depth=" << depth 
00662            << " subdet=" << sdName.Data() << " detId=" << detId << endl << endl; 
00663       return kFALSE;    
00664     }
00665 
00666     phiSymCor[HcalDetId(sd, iEta, iPhi, depth)] = value;
00667   }
00668 
00669   return kTRUE;
00670 }
00671 
00672 
00673 void hcalCalib::makeTextFile() {
00674 
00675   //{ HcalEmpty=0, HcalBarrel=1, HcalEndcap=2, HcalOuter=3, HcalForward=4, HcalTriggerTower=5, HcalOther=7 };
00676   
00677   TString sdName[8] = {"EMPTY", "HB", "HE", "HO", "HF", "TRITWR", "UNDEF", "OTHER"}; 
00678 
00679   FILE *constFile = fopen(OUTPUT_COR_COEF_FILENAME.Data(), "w+"); 
00680 
00681   // header of the constants file
00682   fprintf(constFile, "%1s%16s%16s%16s%16s%9s%11s\n", 
00683           "#", "eta", "phi", "depth", "det", "value", "DetId");
00684 
00685   // Order loops to produce sequence of constants as for phi symmetry 
00686   
00687   for(Int_t sd=1; sd<=4; sd++) {
00688 
00689     for(Int_t e=1; e<=50; e++) {
00690       Int_t eta;
00691 
00692       for(Int_t side=0; side<2; side++) {
00693         eta = (side==0)? -e : e; //ta *= (-1);
00694 
00695         for(Int_t phi=1; phi<=72; phi++) {
00696 
00697           for(Int_t d=1; d<5; d++) {
00698 
00699             if (!HcalDetId::validDetId(HcalSubdetector(sd), eta, phi, d)) continue;
00700             HcalDetId id(HcalSubdetector(sd), eta, phi, d);
00701             Float_t corrFactor = 1.0;
00702             
00703 
00704             if (abs(eta)>=CALIB_ABS_IETA_MIN && abs(eta)<=CALIB_ABS_IETA_MAX && HcalSubdetector(sd)!=HcalOuter) {
00705             //      if (abs(eta)>=CALIB_ABS_IETA_MIN && abs(eta)<=22 && HcalSubdetector(sd)!=HcalOuter) {
00706 
00707 
00708               // need some care when depths were summed for iEta=16 =>
00709               // the coeficients are saved in depth 1 of HB: affects
00710               Int_t subdetInd = sd;
00711               
00712               if (abs(eta)==16 && HcalSubdetector(sd)==HcalEndcap && SUM_DEPTHS) {
00713                 subdetInd = 1;
00714               }
00715 
00716               if (CALIB_METHOD=="L3" || CALIB_METHOD=="L3_AND_MTRX_INV") {
00717 
00718 
00719                 if (SUM_DEPTHS && COMBINE_PHI) corrFactor = solution[HcalDetId(HcalSubdetector(subdetInd), eta, 1, 1)];
00720                 else if (SUM_DEPTHS) corrFactor = solution[HcalDetId(HcalSubdetector(subdetInd), eta, phi, 1)];
00721                 else if (COMBINE_PHI) corrFactor = solution[HcalDetId(HcalSubdetector(sd), eta, 1, d)];
00722                 else corrFactor = solution[HcalDetId(HcalSubdetector(sd), eta, phi, d)];
00723 
00724 
00725                 // Remark: a new case was added (sumSmallDepths) where the first two depths in towers 15,16 
00726                 // are summed and stored in  depth 1.
00727                 // For now we create the correction coef for depth 2 (set equal to depth 1)
00728                 // after the call to the L3 minimizer so that this case is also handled without modifying the
00729                 // logic above. Probably it is better to move it here? 
00730                 
00731 
00732               }// L3
00733 
00734               else if (CALIB_METHOD=="MATRIX_INV_OF_ETA_AVE") {
00735                 corrFactor = iEtaCoefMap[eta];
00736               }
00737 
00738               else if (CALIB_METHOD=="ISO_TRK_PHI_SYM") {
00739                 corrFactor = solution[HcalDetId(HcalSubdetector(sd), eta, phi, d)];
00740               }
00741 
00742 
00743             } // if within the calibration range
00744             
00745             fprintf(constFile, "%17i%16i%16i%16s%9.5f%11X\n", 
00746                     eta, phi, d, sdName[sd].Data(), corrFactor, id.rawId());
00747 
00748           }
00749         }
00750       }
00751     }
00752   }
00753    
00754   return;
00755 }
00756