File indexing completed on 2025-12-19 09:18:41
0001 #include "LiteCaloEval.h"
0002
0003 #include <calobase/RawTower.h>
0004 #include <calobase/RawTowerContainer.h>
0005 #include <calobase/RawTowerGeomContainer.h>
0006 #include <calobase/TowerInfo.h>
0007 #include <calobase/TowerInfoContainer.h>
0008
0009 #include <calotrigger/TriggerAnalyzer.h>
0010
0011 #include <fun4all/Fun4AllReturnCodes.h>
0012 #include <fun4all/SubsysReco.h>
0013
0014 #include <phool/getClass.h>
0015 #include <phool/phool.h>
0016
0017 #include <RtypesCore.h> // for Double_t
0018 #include <TCanvas.h>
0019 #include <TF1.h>
0020 #include <TFile.h>
0021 #include <TGraph.h>
0022 #include <TGraphErrors.h>
0023 #include <TH1.h>
0024 #include <TH2.h>
0025 #include <TH3.h>
0026 #include <TLatex.h>
0027 #include <TLegend.h>
0028 #include <TStyle.h>
0029 #include <TSystem.h>
0030
0031 #include <cmath>
0032 #include <cstdlib>
0033 #include <format>
0034 #include <iostream>
0035 #include <map> // for _Rb_tree_const_iterator
0036 #include <utility> // for pair
0037
0038 class RawTowerGeom;
0039
0040 namespace
0041 {
0042 TGraph *LCE_grff{nullptr};
0043
0044
0045 double LCE_fitf(const Double_t *x, const Double_t *par)
0046 {
0047 return par[0] * LCE_grff->Eval(x[0] * par[1], nullptr, "S");
0048 }
0049 }
0050
0051
0052 LiteCaloEval::LiteCaloEval(const std::string &name, const std::string &caloname, const std::string &filename)
0053 : SubsysReco(name)
0054 , _caloname(caloname)
0055 , _filename(filename)
0056 {
0057 }
0058
0059
0060 int LiteCaloEval::InitRun(PHCompositeNode * )
0061 {
0062 std::cout << "In InitRun " << std::endl;
0063
0064
0065 if (calotype == LiteCaloEval::NONE)
0066 {
0067 std::cout << Name() << ": No Calo Type set" << std::endl;
0068 gSystem->Exit(1);
0069 }
0070
0071 _ievent = 0;
0072
0073 trigAna = new TriggerAnalyzer();
0074
0075 cal_output = new TFile(_filename.c_str(), "RECREATE");
0076
0077 h_event = new TH1F("h_event", "", 1, 0, 1);
0078
0079 if (calotype == LiteCaloEval::HCALIN)
0080 {
0081 hcalin_energy_eta = new TH2F("hcalin_energy_eta", "hcalin energy eta", 100, 0, 10, 24, -0.5, 23.5);
0082 hcalin_e_eta_phi = new TH3F("hcalin_e_eta_phi", "hcalin e eta phi", 60, 0, 6, 24, -0.5, 23.5, 64, -0.5, 63.5);
0083
0084
0085 for (int i = 0; i < 24; i++)
0086 {
0087 for (int j = 0; j < 64; j++)
0088 {
0089 std::string hist_name = "hcal_in_eta_" + std::to_string(i) + "_phi_" + std::to_string(j);
0090
0091 hcal_in_eta_phi[i][j] = new TH1F(hist_name.c_str(), "Hcal_in_energy", 40000, 0, 4);
0092 hcal_in_eta_phi[i][j]->SetXTitle("Energy [GeV]");
0093 }
0094 }
0095
0096
0097 for (int i = 0; i < 25; i++)
0098 {
0099 std::string hist_name = "hcalin_eta_" + std::to_string(i);
0100
0101 if (i < 24)
0102 {
0103 hcalin_eta[i] = new TH1F(hist_name.c_str(), "hcalin eta's", 4000, 0, 4.);
0104 hcalin_eta[i]->SetXTitle("Energy [GeV]");
0105 }
0106 else
0107 {
0108 hcalin_eta[i] = new TH1F(hist_name.c_str(), "hcalin eta's", 40000, 0, 4.);
0109 hcalin_eta[i]->SetXTitle("Energy [GeV]");
0110 }
0111 }
0112 }
0113
0114 else if (calotype == LiteCaloEval::HCALOUT)
0115 {
0116 hcalout_energy_eta = new TH2F("hcalout_energy_eta", "hcalout energy eta", 100, 0, 10, 24, 0.5, 23.5);
0117 hcalout_e_eta_phi = new TH3F("hcalout_e_eta_phi", "hcalout e eta phi", 100, 0, 10, 24, -0.5, 23.5, 64, -0.5, 63.5);
0118
0119
0120 for (int i = 0; i < 24; i++)
0121 {
0122 for (int j = 0; j < 64; j++)
0123 {
0124 std::string hist_name = "hcal_out_eta_" + std::to_string(i) + "_phi_" + std::to_string(j);
0125
0126 hcal_out_eta_phi[i][j] = new TH1F(hist_name.c_str(), "Hcal_out energy", 10000, 0, 10);
0127 hcal_out_eta_phi[i][j]->SetXTitle("Energy [GeV]");
0128 }
0129 }
0130
0131
0132 for (int i = 0; i < 25; i++)
0133 {
0134 std::string hist_name = "hcalout_eta_" + std::to_string(i);
0135 if (i < 24)
0136 {
0137 hcalout_eta[i] = new TH1F(hist_name.c_str(), "hcalout eta's", 10000, 0, 10);
0138 hcalout_eta[i]->SetXTitle("Energy [GeV]");
0139 }
0140 else
0141 {
0142 hcalout_eta[i] = new TH1F(hist_name.c_str(), "hcalout eta's", 100000, 0, 10);
0143 hcalout_eta[i]->SetXTitle("Energy [GeV]");
0144 }
0145 }
0146 }
0147
0148 else if (calotype == LiteCaloEval::CEMC)
0149 {
0150
0151 for (int i = 0; i < 96; i++)
0152 {
0153 for (int j = 0; j < 256; j++)
0154 {
0155 std::string hist_name = "emc_ieta" + std::to_string(i) + "_phi" + std::to_string(j);
0156
0157 cemc_hist_eta_phi[i][j] = new TH1F(hist_name.c_str(), "Hist_ieta_phi_leaf(e)", 400, 0, 2);
0158 cemc_hist_eta_phi[i][j]->SetXTitle("Energy [GeV]");
0159 }
0160 }
0161
0162
0163 for (int i = 0; i < 97; i++)
0164 {
0165 gStyle->SetOptFit(1);
0166 std::string b = "eta_" + std::to_string(i);
0167 eta_hist[i] = new TH1F(b.c_str(), "eta and all phi's", 400, 0, 2);
0168 eta_hist[i]->SetXTitle("Energy [GeV]");
0169 }
0170
0171
0172 energy_eta_hist = new TH2F("energy_eta_hist", "energy eta and all phi", 100, 0, 10, 96, -0.5, 95.5);
0173
0174
0175 e_eta_phi = new TH3F("e_eta_phi", "e v eta v phi", 100, 0, 10, 96, -0.5, 95.5, 256, -0.5, 255.5);
0176 }
0177
0178 return Fun4AllReturnCodes::EVENT_OK;
0179 }
0180
0181
0182 int LiteCaloEval::process_event(PHCompositeNode *topNode)
0183 {
0184 if (_ievent % 100 == 0)
0185 {
0186 std::cout << "LiteCaloEval::process_event(PHCompositeNode *topNode) Processing Event " << _ievent << std::endl;
0187 }
0188
0189
0190 trigAna->decodeTriggers(topNode);
0191
0192 if (reqMinBias && trigAna->didTriggerFire(12) == false)
0193 {
0194 _ievent++;
0195 return Fun4AllReturnCodes::EVENT_OK;
0196 }
0197
0198 h_event->Fill(0);
0199
0200
0201
0202 std::string towernode = "TOWER_CALIB_" + _caloname;
0203 RawTowerContainer *towers = nullptr;
0204 RawTowerGeomContainer *towergeom = nullptr;
0205
0206
0207 RawTowerContainer::ConstRange begin_end;
0208 RawTowerContainer::ConstIterator rtiter;
0209
0210 if (m_UseTowerInfo < 1)
0211 {
0212 towers = findNode::getClass<RawTowerContainer>(topNode, towernode);
0213 if (!towers)
0214 {
0215 std::cout << PHWHERE << " ERROR: Can't find " << towernode << std::endl;
0216 exit(-1);
0217 }
0218
0219
0220 std::string towergeomnode = "TOWERGEOM_" + _caloname;
0221 towergeom = findNode::getClass<RawTowerGeomContainer>(topNode, towergeomnode);
0222 if (!towergeom)
0223 {
0224 std::cout << PHWHERE << " ERROR: Can't find " << towergeomnode << std::endl;
0225 exit(-1);
0226 }
0227
0228 begin_end = towers->getTowers();
0229 rtiter = begin_end.first;
0230 }
0231
0232 if (mode && _ievent < 2)
0233 {
0234 std::cout << "mode is set " << std::endl;
0235 }
0236
0237 TowerInfoContainer *towerinfos = nullptr;
0238
0239
0240 if (m_UseTowerInfo)
0241 {
0242 towernode = _inputnodename;
0243
0244 towerinfos = findNode::getClass<TowerInfoContainer>(topNode, towernode);
0245
0246 if (!towerinfos)
0247 {
0248 std::cout << PHWHERE << " ERROR: Can't find " << towernode << std::endl;
0249 exit(-1);
0250 }
0251 }
0252
0253
0254 unsigned int nchannels = -1;
0255
0256 if (m_UseTowerInfo)
0257 {
0258 nchannels = towerinfos->size();
0259 }
0260 else
0261 {
0262 nchannels = towers->size();
0263 }
0264
0265 TowerInfo *tower_info;
0266 RawTower *tower;
0267 RawTowerGeom *tower_geom;
0268
0269 for (unsigned int channel = 0; channel < nchannels; channel++)
0270 {
0271 if (!m_UseTowerInfo && rtiter == begin_end.second)
0272 {
0273 std::cout << "ERROR: Recheck iteration process with rtiter variable" << std::endl;
0274 }
0275
0276 float e = -1.0;
0277 unsigned int ieta = 99999;
0278 unsigned int iphi = 99999;
0279
0280 if (m_UseTowerInfo)
0281 {
0282 tower_info = towerinfos->get_tower_at_channel(channel);
0283 unsigned int towerkey = towerinfos->encode_key(channel);
0284
0285 e = tower_info->get_energy();
0286 ieta = towerinfos->getTowerEtaBin(towerkey);
0287 iphi = towerinfos->getTowerPhiBin(towerkey);
0288
0289 if (!tower_info->get_isGood())
0290 {
0291 continue;
0292 }
0293 }
0294
0295 else
0296 {
0297 tower = rtiter->second;
0298
0299 tower_geom = towergeom->get_tower_geometry(tower->get_id());
0300
0301 if (!tower_geom)
0302 {
0303 std::cout << PHWHERE << " ERROR: Can't find tower geometry for this tower hit: ";
0304 tower->identify();
0305 exit(-1);
0306 }
0307
0308 e = tower->get_energy();
0309 ieta = tower->get_bineta();
0310 iphi = tower->get_binphi();
0311
0312 }
0313
0314 if (ieta > 95 || iphi > 256)
0315 {
0316
0317 std::cout << "ieta or iphi not set correctly/ was less than 0 " << std::endl;
0318 break;
0319 }
0320
0321 if (e <= 0)
0322 {
0323 continue;
0324 }
0325
0326 if (calotype == LiteCaloEval::CEMC)
0327 {
0328
0329 if (mode)
0330 {
0331 int ket = ieta / 8;
0332 int llet = ket % 6;
0333 int pket = iphi / 8;
0334 int ppkket = pket % 3;
0335
0336 e *= 0.88 + llet * 0.04 - 0.01 + 0.01 * ppkket;
0337 }
0338
0339 cemc_hist_eta_phi[ieta][iphi]->Fill(e);
0340
0341 eta_hist[96]->Fill(e);
0342
0343 eta_hist[ieta]->Fill(e);
0344
0345 energy_eta_hist->Fill(e, ieta);
0346
0347 e_eta_phi->Fill(e, ieta, iphi);
0348 }
0349
0350 else if (calotype == LiteCaloEval::HCALOUT)
0351 {
0352 if (mode)
0353 {
0354 int ket = ieta / 2;
0355
0356 int llet = ket % 6;
0357 e *= 0.945 + llet * 0.02;
0358 int pket = iphi / 4;
0359 if (pket % 2 == 0)
0360 {
0361 e *= 1.03;
0362 }
0363 }
0364
0365 hcal_out_eta_phi[ieta][iphi]->Fill(e);
0366
0367 hcalout_eta[24]->Fill(e);
0368
0369 hcalout_eta[ieta]->Fill(e);
0370
0371 hcalout_energy_eta->Fill(e, ieta);
0372
0373 hcalout_e_eta_phi->Fill(e, ieta, iphi);
0374 }
0375
0376 else if (calotype == LiteCaloEval::HCALIN)
0377 {
0378 if (mode)
0379 {
0380 int ket = ieta / 2;
0381
0382 int llet = ket % 6;
0383 e *= 0.945 + llet * 0.02;
0384 int pket = iphi / 4;
0385 if (pket % 2 == 0)
0386 {
0387 e *= 1.03;
0388 }
0389 }
0390
0391 hcal_in_eta_phi[ieta][iphi]->Fill(e);
0392
0393 hcalin_eta[24]->Fill(e);
0394
0395 hcalin_eta[ieta]->Fill(e);
0396
0397 hcalin_energy_eta->Fill(e, ieta);
0398
0399 hcalin_e_eta_phi->Fill(e, ieta, iphi);
0400 }
0401
0402 if (!m_UseTowerInfo)
0403 {
0404 ++rtiter;
0405 }
0406
0407 }
0408
0409 _ievent++;
0410
0411 return Fun4AllReturnCodes::EVENT_OK;
0412 }
0413
0414
0415 int LiteCaloEval::End(PHCompositeNode * )
0416 {
0417 cal_output->cd();
0418
0419 std::cout << " writing lite calo file" << std::endl;
0420
0421 cal_output->Write();
0422
0423 return Fun4AllReturnCodes::EVENT_OK;
0424 }
0425
0426
0427 void LiteCaloEval::Get_Histos(const std::string &infile, const std::string &outfile)
0428 {
0429 std::cout << "Getting histograms... " << std::endl;
0430
0431 if (infile.empty())
0432 {
0433 std::cout << "Need an input file to run LiteCaloEval::Get_Histos" << std::endl;
0434 exit(0);
0435 }
0436
0437 if (!outfile.empty())
0438 {
0439 std::string ts = "cp ";
0440 ts += infile;
0441 ts += " ";
0442 ts += outfile;
0443 gSystem->Exec(ts.c_str());
0444 f_temp = new TFile(outfile.c_str(), "UPDATE");
0445 }
0446
0447 else
0448 {
0449 f_temp = new TFile(infile.c_str());
0450 }
0451
0452 int max_ieta = 96;
0453 if (calotype != LiteCaloEval::CEMC)
0454 {
0455 max_ieta = 24;
0456 }
0457
0458 int max_iphi = 256;
0459 if (calotype != LiteCaloEval::CEMC)
0460 {
0461 max_iphi = 64;
0462 }
0463
0464
0465 for (int i = 0; i < max_ieta + 1; i++)
0466 {
0467 std::string b = "eta_" + std::to_string(i);
0468
0469 if (calotype == LiteCaloEval::HCALOUT)
0470 {
0471 b.insert(0, "hcalout_");
0472 }
0473 else if (calotype == LiteCaloEval::HCALIN)
0474 {
0475 b.insert(0, "hcalin_");
0476 }
0477
0478
0479 TH1 *heta_temp{nullptr};
0480 f_temp->GetObject(b.c_str(), heta_temp);
0481
0482 if (i == 0)
0483 {
0484 std::cout << "Old bin width for eta slice " << heta_temp->GetBinWidth(1) << std::endl;
0485 }
0486
0487
0488 heta_temp->Rebin(binwidth / heta_temp->GetBinWidth(1));
0489
0490 if (i == 0)
0491 {
0492 std::cout << "New bin width for eta slice " << heta_temp->GetBinWidth(1) << std::endl;
0493 }
0494
0495 if (!heta_temp && i == 0)
0496 {
0497 std::cout << " warning hist " << b << " not found" << std::endl;
0498 }
0499
0500
0501 eta_hist[i] = heta_temp;
0502
0503 if (calotype == LiteCaloEval::HCALOUT)
0504 {
0505 hcalout_eta[i] = heta_temp;
0506 }
0507 else if (calotype == LiteCaloEval::HCALIN)
0508 {
0509 hcalin_eta[i] = heta_temp;
0510 }
0511
0512 if (!(i < max_ieta))
0513 {
0514 continue;
0515 }
0516
0517
0518 for (int j = 0; j < max_iphi; j++)
0519 {
0520
0521 std::string hist_name_p = "emc_ieta" + std::to_string(i) + "_phi" + std::to_string(j);
0522
0523 if (calotype == LiteCaloEval::HCALOUT)
0524 {
0525 hist_name_p = "hcal_out_eta_" + std::to_string(i) + "_phi_" + std::to_string(j);
0526 }
0527 else if (calotype == LiteCaloEval::HCALIN)
0528 {
0529 hist_name_p = "hcal_in_eta_" + std::to_string(i) + "_phi_" + std::to_string(j);
0530 }
0531
0532
0533 TH1 *heta_tempp{nullptr};
0534 f_temp->GetObject(hist_name_p.c_str(), heta_tempp);
0535
0536 if (i == 0 && j == 0)
0537 {
0538 std::cout << "Old bin width for towers " << heta_tempp->GetBinWidth(1) << std::endl;
0539 }
0540
0541
0542 heta_tempp->Rebin(binwidth / heta_tempp->GetBinWidth(1));
0543
0544 if (i == 0 && j == 0)
0545 {
0546 std::cout << "New bin width for towers " << heta_tempp->GetBinWidth(1) << std::endl;
0547 }
0548
0549 if (!heta_tempp && i == 0)
0550 {
0551 std::cout << " warning hist " << hist_name_p.c_str() << " not found" << std::endl;
0552 }
0553
0554
0555 cemc_hist_eta_phi[i][j] = heta_tempp;
0556
0557 if (calotype == LiteCaloEval::HCALOUT)
0558 {
0559 hcal_out_eta_phi[i][j] = heta_tempp;
0560 }
0561 else if (calotype == LiteCaloEval::HCALIN)
0562 {
0563 hcal_in_eta_phi[i][j] = heta_tempp;
0564 }
0565 }
0566 }
0567
0568
0569
0570
0571
0572
0573
0574 std::cout << "Grabbed all histograms." << std::endl;
0575
0576 }
0577
0578 void LiteCaloEval::FitRelativeShifts(LiteCaloEval *ref_lce, int modeFitShifts)
0579 {
0580 bool onlyEta = false;
0581
0582 if (fitmin < 0.001)
0583 {
0584 fitmin = 0.15;
0585 }
0586 if (fitmax < 0.001)
0587 {
0588 fitmax = 1.3;
0589 }
0590
0591 float par_value[96] = {0};
0592 float par_err[96] = {0};
0593 float eta_value[96] = {0};
0594 float eta_err[96] = {0};
0595
0596 if (f_temp)
0597 {
0598 f_temp->cd();
0599 }
0600
0601
0602 TF1 *f1 = new TF1("myexpo", LCE_fitf, 0.1, 10, 2);
0603 f1->SetParameters(1.0, 1.0);
0604
0605
0606 if (modeFitShifts % 10 == 1)
0607 {
0608 onlyEta = true;
0609 std::cout << "Running only eta slice mode....." << std::endl;
0610 }
0611
0612
0613 int nsmooth = 1;
0614
0615
0616 int addSmooth = modeFitShifts % 100 / 10;
0617 if (addSmooth)
0618 {
0619 nsmooth += addSmooth;
0620 }
0621
0622
0623 bool flag_fit_rings = false;
0624
0625
0626 if (modeFitShifts % 1000 / 100 == 1)
0627 {
0628 flag_fit_rings = true;
0629 }
0630
0631 int max_ieta = 96;
0632 if (calotype != LiteCaloEval::CEMC)
0633 {
0634 max_ieta = 24;
0635 }
0636
0637 int max_iphi = 256;
0638 if (calotype != LiteCaloEval::CEMC)
0639 {
0640 max_iphi = 64;
0641 }
0642
0643
0644 TH2 *corrPat = new TH2F("corrPat", "", max_ieta, 0, max_ieta, max_iphi, 0, max_iphi);
0645 corrPat->SetXTitle("#eta bin");
0646 corrPat->SetYTitle("#phi bin");
0647
0648 TH2 *h2_failQA = new TH2F("h2_failQA", "", max_ieta, 0, max_ieta, max_iphi, 0, max_iphi);
0649
0650
0651 TH1 *gainvals = new TH1F("gainvals", "Towerslope Correction Values", 10000, 0, 10);
0652 gainvals->SetXTitle("Gain Shift Value");
0653 gainvals->SetYTitle("Counts");
0654
0655
0656 TH1 *h_gainErr = new TH1F("h_gainErr", "Towerslope Corrections Errors", 1000, 0, 1);
0657 h_gainErr->SetXTitle("error");
0658 h_gainErr->SetYTitle("Counts");
0659
0660 int minbin = 0;
0661 int maxbin = max_ieta;
0662
0663 if (m_myminbin > -1)
0664 {
0665 minbin = m_myminbin;
0666 }
0667 if (m_mymaxbin > -1)
0668 {
0669 maxbin = m_mymaxbin;
0670 }
0671
0672
0673 TH1 *hnewf = nullptr;
0674
0675
0676 for (int i = minbin; i < maxbin; i++)
0677 {
0678 std::cout << "===============================" << std::endl;
0679 std::cout << "Fitting towers in eta slice " << i << std::endl;
0680
0681
0682 double ieta_gain;
0683 double ieta_gain_err;
0684
0685
0686 std::string myClnm = "newhc_eta" + std::to_string(i);
0687
0688
0689 int iik = i;
0690
0691
0692 TH1 *cleanEtaRef = nullptr;
0693 std::string cleanEta = "cleanEtaRef_";
0694
0695 if (calotype == LiteCaloEval::CEMC)
0696 {
0697 if (flag_fit_rings == true)
0698 {
0699 cleanEta += std::to_string(iik);
0700
0701 cleanEtaRef = (TH1 *) eta_hist[iik]->Clone(cleanEta.c_str());
0702 }
0703
0704 else
0705 {
0706 hnewf = (TH1 *) ref_lce->eta_hist[iik]->Clone(myClnm.c_str());
0707 }
0708 }
0709
0710 else if (calotype == LiteCaloEval::HCALOUT)
0711 {
0712 if (flag_fit_rings == true)
0713 {
0714 cleanEta += std::to_string(iik);
0715
0716 cleanEtaRef = (TH1 *) hcalout_eta[iik]->Clone(cleanEta.c_str());
0717
0718
0719 if (i < 4 || i > 19)
0720 {
0721 for (int phiCH = 14; phiCH < 20; phiCH++)
0722 {
0723 cleanEtaRef->Add(hcal_out_eta_phi[i][phiCH], -1.0);
0724 }
0725 }
0726 }
0727
0728 else
0729 {
0730 hnewf = (TH1 *) ref_lce->hcalout_eta[iik]->Clone(myClnm.c_str());
0731 }
0732 }
0733
0734 else if (calotype == LiteCaloEval::HCALIN)
0735 {
0736 if (flag_fit_rings == true)
0737 {
0738 cleanEta += std::to_string(iik);
0739
0740 cleanEtaRef = (TH1 *) hcalin_eta[iik]->Clone(cleanEta.c_str());
0741 }
0742
0743 else
0744 {
0745 hnewf = (TH1 *) ref_lce->hcalin_eta[iik]->Clone(myClnm.c_str());
0746 }
0747 }
0748
0749 if (flag_fit_rings == true)
0750 {
0751 cleanEtaRef->Smooth(nsmooth);
0752
0753 LCE_grff = new TGraph(cleanEtaRef);
0754 }
0755
0756 else
0757 {
0758 hnewf->Smooth(nsmooth);
0759
0760 LCE_grff = new TGraph(hnewf);
0761 }
0762
0763
0764 TF1 *f2f = nullptr;
0765
0766
0767 if (calotype == LiteCaloEval::CEMC)
0768 {
0769 if (nsmooth > 1)
0770 {
0771 eta_hist[i]->Smooth(nsmooth);
0772 }
0773
0774 eta_hist[i]->Fit("myexpo", "Q", "", fitmin, fitmax);
0775
0776 f2f = eta_hist[i]->GetFunction("myexpo");
0777 }
0778
0779 else if (calotype == LiteCaloEval::HCALOUT)
0780 {
0781 if (nsmooth > 1)
0782 {
0783 hcalout_eta[i]->Smooth(nsmooth);
0784 }
0785
0786 hcalout_eta[i]->Fit("myexpo", "Q", "", fitmin, fitmax);
0787
0788 f2f = hcalout_eta[i]->GetFunction("myexpo");
0789 }
0790
0791 else if (calotype == LiteCaloEval::HCALIN)
0792 {
0793 if (nsmooth > 1)
0794 {
0795 hcalin_eta[i]->Smooth(nsmooth);
0796 }
0797
0798 hcalin_eta[i]->Fit("myexpo", "Q", "", fitmin, fitmax);
0799
0800 f2f = hcalin_eta[i]->GetFunction("myexpo");
0801 }
0802
0803 if (!f2f)
0804 {
0805 std::cout << "Warning, f2f is null!" << std::endl;
0806 exit(-1);
0807 }
0808 else
0809 {
0810 ieta_gain = f2f->GetParameter(1);
0811 }
0812
0813 ieta_gain_err = f2f->GetParError(1);
0814
0815 par_value[i] = ieta_gain;
0816 par_err[i] = ieta_gain_err;
0817 eta_value[i] = i;
0818 eta_err[i] = 0.0;
0819
0820
0821 if (onlyEta == true)
0822 {
0823 continue;
0824 }
0825
0826
0827
0828 if (doQA)
0829 {
0830 if (flag_fit_rings == true)
0831 {
0832 if (calotype == LiteCaloEval::CEMC)
0833 {
0834 for (int j = 0; j < max_iphi; j++)
0835 {
0836 cleanEtaRef->Add(cemc_hist_eta_phi[i][j], -1.0);
0837
0838 bool qa_res = spec_QA(cemc_hist_eta_phi[i][j], cleanEtaRef, true);
0839 if (qa_res == false)
0840 {
0841 std::cout << "failed QA " << i << "," << j << std::endl;
0842 cemc_hist_eta_phi[i][j]->Reset();
0843 h2_failQA->Fill(i, j);
0844 }
0845 else
0846 {
0847 cleanEtaRef->Add(cemc_hist_eta_phi[i][j], 1.0);
0848 }
0849 }
0850 }
0851 }
0852 }
0853
0854
0855
0856
0857 for (int j = 0; j < max_iphi; j++)
0858 {
0859
0860 std::string myClnmp = "newhc_eta" + std::to_string((1000 * (i + 2)) + j);
0861
0862
0863 TH1 *hnewfp = nullptr;
0864
0865
0866
0867 bool _isChimney = chk_isChimney(i, j);
0868
0869 if (calotype == LiteCaloEval::CEMC)
0870 {
0871
0872 if (!(cemc_hist_eta_phi[i][j]->GetEntries()))
0873 {
0874 std::cout << "No entries in EMCAL tower histogram (" << i << "," << j << "). Skipping fitting." << std::endl;
0875 continue;
0876 }
0877
0878 if (flag_fit_rings == true)
0879 {
0880 cleanEtaRef->Add(cemc_hist_eta_phi[i][j], -1.0);
0881 }
0882 else
0883 {
0884 hnewfp = (TH1 *) ref_lce->cemc_hist_eta_phi[i][j]->Clone(myClnmp.c_str());
0885 }
0886 }
0887
0888 else if (calotype == LiteCaloEval::HCALOUT)
0889 {
0890
0891 if (!(hcal_out_eta_phi[i][j]->GetEntries()))
0892 {
0893 std::cout << "No entries in OHCAL tower histogram (" << i << "," << j << "). Skipping fitting." << std::endl;
0894 continue;
0895 }
0896
0897 if (flag_fit_rings == true)
0898 {
0899
0900
0901 if (!_isChimney)
0902 {
0903 cleanEtaRef->Add(hcal_out_eta_phi[i][j], -1.0);
0904 }
0905 }
0906
0907 else
0908 {
0909 hnewfp = (TH1 *) ref_lce->hcal_out_eta_phi[i][j]->Clone(myClnmp.c_str());
0910 }
0911 }
0912
0913 else if (calotype == LiteCaloEval::HCALIN)
0914 {
0915
0916 if (!(hcal_in_eta_phi[i][j]->GetEntries()))
0917 {
0918 std::cout << "No entries in IHCAL tower histogram(" << i << "," << j << "). Skipping fitting." << std::endl;
0919 continue;
0920 }
0921
0922 if (flag_fit_rings == true)
0923 {
0924 cleanEtaRef->Add(hcal_in_eta_phi[i][j], -1.0);
0925 }
0926 else
0927 {
0928 hnewfp = (TH1 *) ref_lce->hcal_in_eta_phi[i][j]->Clone(myClnmp.c_str());
0929 }
0930 }
0931
0932
0933
0934
0935
0936
0937
0938 if (flag_fit_rings == false)
0939 {
0940 hnewfp->Smooth(nsmooth);
0941 LCE_grff = new TGraph(hnewfp);
0942 }
0943
0944 if (flag_fit_rings == true)
0945 {
0946
0947 LCE_grff = new TGraph(cleanEtaRef);
0948 }
0949
0950
0951 TF1 *f2f2 = nullptr;
0952
0953
0954 if (calotype == LiteCaloEval::CEMC)
0955 {
0956 double scaleP0 = cemc_hist_eta_phi[i][j]->Integral(cemc_hist_eta_phi[i][j]->FindBin(fitmin), cemc_hist_eta_phi[i][j]->FindBin(fitmax));
0957
0958 if (flag_fit_rings == 1)
0959 {
0960 scaleP0 /= cleanEtaRef->Integral(cleanEtaRef->FindBin(fitmin), cleanEtaRef->FindBin(fitmax));
0961 }
0962 else
0963 {
0964 scaleP0 /= hnewfp->Integral(hnewfp->FindBin(fitmin), hnewfp->FindBin(fitmax));
0965 }
0966
0967 f1->SetParameters(scaleP0, 1.0);
0968
0969 cemc_hist_eta_phi[i][j]->Fit("myexpo", "Q", "", fitmin, fitmax);
0970
0971 f2f2 = cemc_hist_eta_phi[i][j]->GetFunction("myexpo");
0972
0973
0974 if (flag_fit_rings == true)
0975 {
0976 cleanEtaRef->Add(cemc_hist_eta_phi[i][j], 1.0);
0977 }
0978 }
0979
0980 else if (calotype == LiteCaloEval::HCALOUT)
0981 {
0982 double scaleP0 = hcal_out_eta_phi[i][j]->Integral(hcal_out_eta_phi[i][j]->FindBin(fitmin), hcal_out_eta_phi[i][j]->FindBin(fitmax));
0983
0984 if (flag_fit_rings == 1)
0985 {
0986 scaleP0 /= cleanEtaRef->Integral(cleanEtaRef->FindBin(fitmin), cleanEtaRef->FindBin(fitmax));
0987 }
0988 else
0989 {
0990 scaleP0 /= hnewfp->Integral(hnewfp->FindBin(fitmin), hnewfp->FindBin(fitmax));
0991 }
0992
0993 f1->SetParameters(scaleP0, 1.0);
0994
0995 hcal_out_eta_phi[i][j]->Fit("myexpo", "Q", "", fitmin, fitmax);
0996
0997 f2f2 = hcal_out_eta_phi[i][j]->GetFunction("myexpo");
0998
0999 if (flag_fit_rings == true)
1000 {
1001 if (!_isChimney)
1002 {
1003 cleanEtaRef->Add(hcal_out_eta_phi[i][j], 1.0);
1004 }
1005 }
1006 }
1007
1008 else if (calotype == LiteCaloEval::HCALIN)
1009 {
1010 double scaleP0 = hcal_in_eta_phi[i][j]->Integral(hcal_in_eta_phi[i][j]->FindBin(fitmin), hcal_in_eta_phi[i][j]->FindBin(fitmax));
1011
1012 if (flag_fit_rings == 1)
1013 {
1014 scaleP0 /= cleanEtaRef->Integral(cleanEtaRef->FindBin(fitmin), cleanEtaRef->FindBin(fitmax));
1015 }
1016
1017 else
1018 {
1019 scaleP0 /= hnewfp->Integral(hnewfp->FindBin(fitmin), hnewfp->FindBin(fitmax));
1020 }
1021
1022 f1->SetParameters(scaleP0, 1.0);
1023
1024 hcal_in_eta_phi[i][j]->Fit("myexpo", "Q", "", fitmin, fitmax);
1025
1026 f2f2 = hcal_in_eta_phi[i][j]->GetFunction("myexpo");
1027
1028 if (flag_fit_rings == true)
1029 {
1030 cleanEtaRef->Add(hcal_in_eta_phi[i][j], 1.0);
1031 }
1032 }
1033
1034 float correction = 1;
1035 float corrErr = 1;
1036
1037 if (f2f2)
1038 {
1039 correction = f2f2->GetParameter(1);
1040 corrErr = f2f2->GetParError(1);
1041 }
1042
1043 double errProp = corrErr / (correction * correction);
1044
1045
1046
1047
1048
1049 corrPat->SetBinContent(i + 1, j + 1, 1 / correction);
1050
1051 corrPat->SetBinError(i + 1, j + 1, errProp);
1052
1053 gainvals->Fill(1.0 / correction);
1054
1055 h_gainErr->Fill(errProp);
1056
1057 }
1058
1059 }
1060
1061
1062 TGraphErrors g1(max_ieta, eta_value, par_value, eta_err, par_err);
1063 g1.SetTitle("fitted shifts; eta; p1");
1064 g1.SetMarkerStyle(20);
1065 g1.Draw("ap");
1066 g1.SetName("Fit1_etaout");
1067 g1.Write();
1068
1069 corrPat->Write();
1070 h2_failQA->Write();
1071 gainvals->Write();
1072 h_gainErr->Write();
1073
1074 if (f_temp)
1075 {
1076 f_temp->Write();
1077 }
1078
1079 f_temp->Close();
1080 }
1081
1082 bool LiteCaloEval::chk_isChimney(int ieta, int iphi)
1083 {
1084 if ((ieta < 4 || ieta > 19) && (iphi > 13 && iphi < 20))
1085 {
1086 return true;
1087 }
1088
1089 return false;
1090 }
1091
1092 float LiteCaloEval::spec_QA(TH1 *h_spec, TH1 *h_ref, bool )
1093 {
1094 float norm = h_spec->Integral(h_spec->FindBin(0.4, 1.1), h_spec->FindBin(1.0));
1095 if (norm == 0)
1096 {
1097 std::cout << "no data for QA" << std::endl;
1098 return false;
1099 }
1100 TF1 *f_exp = new TF1("f_exp", "[0]*TMath::Exp(-1.0*x*[1])", 0.5, 1);
1101 h_spec->Fit(f_exp, "QNR", "", 0.4, 1.);
1102 float spec_p1 = f_exp->GetParameter(1);
1103 h_ref->Fit(f_exp, "QNR", "", 0.4, 1.);
1104 float ref_p1 = f_exp->GetParameter(1);
1105 float calib = ref_p1 / spec_p1;
1106 delete f_exp;
1107 return calib;
1108 }
1109
1110 bool LiteCaloEval::spec_QA(TH1 *h_spec, TH1 *h_ref)
1111 {
1112 float calib = spec_QA(h_spec, h_ref, true);
1113 if (calib > 1.2 || calib < 0.8)
1114 {
1115 return false;
1116 }
1117
1118 return true;
1119 }
1120
1121 void LiteCaloEval::plot_cemc(const std::string &path)
1122 {
1123 TH2 *h_fail = new TH2F("h_fail", "", 96, 0, 96, 256, 0, 256);
1124 TH2 *h_corrPat = new TH2F("corrPat", "", 96, 0, 96, 256, 0, 256);
1125 TH2 *h_hits = new TH2F("h_hits", "", 96, 0, 96, 256, 0, 256);
1126 TH2 *h_expCalib = new TH2F("h_expCalib", "", 96, 0, 96, 256, 0, 256);
1127 h_hits->SetXTitle("#it{#eta}_{i}");
1128 h_hits->SetYTitle("#it{#phi}_{i}");
1129 TH2 *h_hits_clean = new TH2F("h_hits_clean", "", 96, 0, 96, 256, 0, 256);
1130 h_hits_clean->SetXTitle("#it{#eta}_{i}");
1131 h_hits_clean->SetYTitle("#it{#phi}_{i}");
1132
1133 for (int ieta = 0; ieta < 96; ieta++)
1134 {
1135 for (int iphi = 0; iphi < 256; iphi++)
1136 {
1137 h_corrPat->SetBinContent(ieta + 1, iphi + 1, 1);
1138 cemc_hist_eta_phi[ieta][iphi]->Rebin(5);
1139
1140 int binhit = cemc_hist_eta_phi[ieta][iphi]->FindBin(0.25);
1141 int binMax = cemc_hist_eta_phi[ieta][iphi]->GetNbinsX();
1142 float hits = cemc_hist_eta_phi[ieta][iphi]->Integral(binhit, binMax + 1);
1143 h_hits->SetBinContent(ieta + 1, iphi + 1, hits);
1144 h_hits->SetBinError(ieta + 1, iphi + 1, std::sqrt(hits));
1145 h_hits_clean->SetBinContent(ieta + 1, iphi + 1, hits);
1146 h_hits_clean->SetBinError(ieta + 1, iphi + 1, std::sqrt(hits));
1147 int bin = cemc_hist_eta_phi[ieta][iphi]->FindBin(0.2);
1148 if (cemc_hist_eta_phi[ieta][iphi]->GetBinContent(bin) == 0)
1149 {
1150 continue;
1151 }
1152 float calib = spec_QA(cemc_hist_eta_phi[ieta][iphi], eta_hist[ieta], true);
1153 h_expCalib->SetBinContent(ieta + 1, iphi + 1, calib);
1154 if (!spec_QA(cemc_hist_eta_phi[ieta][iphi], eta_hist[ieta], true))
1155 {
1156 h_fail->Fill(ieta, iphi);
1157 h_corrPat->SetBinContent(ieta + 1, iphi + 1, 0);
1158 h_hits_clean->SetBinContent(ieta + 1, iphi + 1, 0);
1159 }
1160 }
1161 }
1162
1163
1164 TH2 *h_hot = new TH2F("h_hot", "", 96, 0, 96, 256, 0, 256);
1165 TH1 *h1_hits[96];
1166 TH1 *h1_hits2[96];
1167 float max = h_hits_clean->GetBinContent(h_hits_clean->GetMaximumBin());
1168 float min = h_hits_clean->GetMinimum();
1169 if (min == 0)
1170 {
1171 min = 1;
1172 }
1173
1174 TCanvas *c3 = new TCanvas("c3", "c3", 1600, 600);
1175
1176 for (int ie = 0; ie < 96; ie++)
1177 {
1178 std::string name1 = "h1_hits" + std::to_string(ie);
1179 std::string name2 = "h1_hits2_" + std::to_string(ie);
1180
1181 h1_hits[ie] = new TH1F(name1.c_str(), "", 200, min, max);
1182 h1_hits2[ie] = new TH1F(name2.c_str(), "", 200, min, max);
1183 for (int iphi = 0; iphi < 256; iphi++)
1184 {
1185 h1_hits[ie]->Fill(h_hits_clean->GetBinContent(ie + 1, iphi + 1));
1186 }
1187 float mean = h1_hits[ie]->GetMean();
1188 float std = h1_hits[ie]->GetStdDev();
1189 for (int iphi = 0; iphi < 256; iphi++)
1190 {
1191 float val = h_hits_clean->GetBinContent(ie + 1, iphi + 1);
1192 if (std::fabs(mean - val) / std > 3)
1193 {
1194 continue;
1195 }
1196 h1_hits2[ie]->Fill(val);
1197 }
1198 mean = h1_hits[ie]->GetMean();
1199 std = h1_hits[ie]->GetStdDev();
1200 for (int iphi = 0; iphi < 256; iphi++)
1201 {
1202 h_hot->SetBinContent(ie + 1, iphi + 1, 0);
1203 float val = h_hits_clean->GetBinContent(ie + 1, iphi + 1);
1204 if (val == 0)
1205 {
1206 continue;
1207 }
1208 if (std::fabs(mean - val) / std > 5)
1209 {
1210 h_hot->SetBinContent(ie + 1, iphi + 1, 1);
1211 h_corrPat->SetBinContent(ie + 1, iphi + 1, 0);
1212 h_hits_clean->SetBinContent(ie + 1, iphi + 1, 0);
1213 }
1214 }
1215 h1_hits[ie]->Write();
1216 h1_hits2[ie]->Write();
1217 if (ie < 90)
1218 {
1219 continue;
1220 }
1221 h1_hits[ie]->Draw("hist same");
1222 h1_hits[ie]->GetYaxis()->SetRangeUser(0.1, 100);
1223 h1_hits[ie]->GetXaxis()->SetTitle("Number of hits");
1224 h1_hits[ie]->GetXaxis()->SetNdivisions(505);
1225 h1_hits2[ie]->Draw("hist same");
1226 h1_hits2[ie]->SetLineColor(kBlue);
1227 }
1228 c3->SaveAs((path + "/hits.pdf").c_str());
1229
1230 TLatex latex;
1231 latex.SetTextFont(42);
1232 latex.SetTextSize(0.04);
1233 latex.SetTextAlign(12);
1234 latex.SetNDC();
1235
1236 int ieta = 0;
1237 for (int iplot = 0; iplot < 10; iplot++)
1238 {
1239 TCanvas *c1 = new TCanvas("c1", "c1", 400, 600);
1240
1241 for (int ie = 0; ie < 10; ie++)
1242 {
1243 if (ieta > 95)
1244 {
1245 continue;
1246 }
1247 for (int iphi = 0; iphi < 256; iphi++)
1248 {
1249
1250 int bin = cemc_hist_eta_phi[ieta][iphi]->FindBin(0.2);
1251 if (cemc_hist_eta_phi[ieta][iphi]->GetBinContent(bin) == 0)
1252 {
1253 continue;
1254 }
1255 cemc_hist_eta_phi[ieta][iphi]->Scale(1.0 / cemc_hist_eta_phi[ieta][iphi]->GetBinContent(bin) * pow(10, ie));
1256
1257 cemc_hist_eta_phi[ieta][iphi]->Draw("same hist");
1258 cemc_hist_eta_phi[ieta][iphi]->SetXTitle("pre-calib tower [GeV]");
1259 cemc_hist_eta_phi[ieta][iphi]->SetYTitle("N_{twr}/N_{twr}(E=0.2 GeV) x10^{ieta}");
1260 cemc_hist_eta_phi[ieta][iphi]->GetXaxis()->SetRangeUser(0, 1);
1261 cemc_hist_eta_phi[ieta][iphi]->GetYaxis()->SetRangeUser(0.01, 1e11);
1262 if (h_corrPat->GetBinContent(ieta + 1, iphi + 1) == 0)
1263 {
1264 cemc_hist_eta_phi[ieta][iphi]->SetLineColor(kBlue);
1265 h_fail->Fill(ieta, iphi);
1266 h_hits_clean->SetBinContent(ieta + 1, iphi + 1, 0);
1267 }
1268 latex.SetTextSize(0.02);
1269 latex.SetTextColor(kBlack);
1270 std::string result = "#eta_{" + std::to_string(ieta) + "}#times10^{" + std::to_string(ie) + "}";
1271 latex.DrawLatex(0.22 + (0.07 * ie), 0.92, result.c_str());
1272 }
1273
1274
1275 int bin = eta_hist[ieta]->FindBin(0.2);
1276 if (eta_hist[ieta]->GetBinContent(bin) == 0)
1277 {
1278 ieta++;
1279 std::cout << "no data for eta=" << ieta << std::endl;
1280 continue;
1281 }
1282 eta_hist[ieta]->Scale(1.0 / eta_hist[ieta]->GetBinContent(bin) * pow(10, ie));
1283 eta_hist[ieta]->Draw("same hist");
1284 eta_hist[ieta]->SetLineColor(kRed);
1285
1286 ieta++;
1287 }
1288
1289 gPad->SetLogy();
1290
1291 latex.SetTextSize(0.04);
1292
1293 latex.SetTextColor(kBlack);
1294 latex.DrawLatex(0.6, 0.87, "#it{#bf{sPHENIX}} Internal");
1295 latex.SetTextColor(kRed);
1296 latex.DrawLatex(0.6, 0.83, "eta integrated");
1297 latex.SetTextColor(kBlue);
1298 latex.DrawLatex(0.6, 0.79, "failed QA");
1299
1300 c1->SaveAs((path + "/tsc_spec" + std::to_string(iplot) + ".pdf").c_str());
1301 delete c1;
1302 }
1303
1304 TCanvas *c2 = new TCanvas("c2", "c2", 600, 600);
1305 h_fail->Draw("COLZ");
1306
1307 latex.SetTextColor(kBlack);
1308 std::string stuff = std::to_string(h_fail->GetEntries()) + " failed QA";
1309 latex.DrawLatex(0.2, 0.87, stuff.c_str());
1310
1311 c2->SaveAs((path + "/failQA.pdf").c_str());
1312
1313 TCanvas *c4 = new TCanvas("c4", "c4", 600, 600);
1314 h_expCalib->Draw("COLZ");
1315 h_expCalib->GetZaxis()->SetRangeUser(0.5, 2);
1316 gPad->SetRightMargin(0.15);
1317
1318 latex.SetTextColor(kBlack);
1319
1320 c4->SaveAs((path + "/expCalib.pdf").c_str());
1321
1322 h_hot->Write();
1323 h_hits->Write();
1324 h_hits_clean->Write();
1325 h_corrPat->Write();
1326 h_expCalib->Write();
1327
1328 f_temp->Close();
1329 }
1330
1331 void LiteCaloEval::draw_spectra(const char *outfile)
1332 {
1333 if (calotype == LiteCaloEval::NONE)
1334 {
1335 std::cout << "Did not enter correct calotype. Exiting macro..." << std::endl;
1336 exit(-1);
1337 }
1338
1339 TFile *fout = new TFile(outfile, "UPDATE");
1340
1341 TH1 *h = nullptr;
1342 std::string histName;
1343
1344 TH1 *h_etaSlice = nullptr;
1345 std::string etaSliceName;
1346
1347 TH1 *h_cln = nullptr;
1348 std::string h_cln_nm;
1349
1350 TH1 *h_es_cln = nullptr;
1351 std::string h_es_nm;
1352
1353 float scale = 1.0;
1354
1355 float fitMin = getFitMin();
1356 float fitMax = getFitMax();
1357
1358 double xaxisRange = 2.0;
1359
1360 if (calotype == LiteCaloEval::HCALIN)
1361 {
1362 xaxisRange = 1.0;
1363 }
1364
1365 if (calotype == LiteCaloEval::HCALIN || calotype == LiteCaloEval::HCALOUT)
1366 {
1367 std::cout << "Drawing hcal tower spectra" << std::endl;
1368
1369 int starteta = 0;
1370 int maxeta = 12;
1371 int power = 23;
1372
1373 TCanvas *c = new TCanvas();
1374 c->Divide(2);
1375 c->SetName("hcal_spectra");
1376
1377 int cntr = 1;
1378
1379
1380 for (int k = 1; k < 3; k++)
1381 {
1382 c->cd(k);
1383 gPad->SetLogy(1);
1384
1385 if (k == 2)
1386 {
1387 power = 47;
1388 starteta = 12;
1389 maxeta = 24;
1390 }
1391
1392 TLegend *t = new TLegend(0.5, 0.8, 0.75, 0.9);
1393 t->SetFillStyle(0);
1394 t->SetBorderSize(0);
1395 t->AddEntry("", std::format("ieta {} - {}", starteta, (maxeta - 1)).c_str());
1396
1397
1398 for (int i = starteta; i < maxeta; i++)
1399 {
1400
1401 if (calotype == LiteCaloEval::HCALOUT)
1402 {
1403 etaSliceName = "hcalout_eta_" + std::to_string(i);
1404 }
1405 else if (calotype == LiteCaloEval::HCALIN)
1406 {
1407 etaSliceName = "hcalin_eta_" + std::to_string(i);
1408 }
1409
1410 fout->GetObject(etaSliceName.c_str(), h_etaSlice);
1411
1412 if (!h_etaSlice)
1413 {
1414 std::cout << "ERROR! Could not get hcal eta slice histogram " << i << " ." << std::endl;
1415 gSystem->Exit(1);
1416 exit(1);
1417 }
1418
1419 if (h_etaSlice->GetEntries() == 0.0)
1420 {
1421 std::cout << "WARNING! hcal eta slice " << i << " has no entries!" << std::endl;
1422 continue;
1423 }
1424
1425 h_es_cln = (TH1 *) h_etaSlice->Clone();
1426
1427 h_es_nm = "cln_es_" + std::to_string(i);
1428
1429 h_es_cln->SetName(h_es_nm.c_str());
1430
1431
1432 for (int j = 0; j < 64; j++)
1433 {
1434 if (calotype == LiteCaloEval::HCALOUT)
1435 {
1436 histName = "hcal_out_eta_" + std::to_string(i) + "_phi_" + std::to_string(j);
1437 }
1438 else if (calotype == LiteCaloEval::HCALIN)
1439 {
1440 histName = "hcal_in_eta_" + std::to_string(i) + "_phi_" + std::to_string(j);
1441 }
1442
1443 fout->GetObject(histName.c_str(), h);
1444
1445 if (!h)
1446 {
1447 std::cout << "ERROR! Could not find tower " << histName << "." << std::endl;
1448 gSystem->Exit(1);
1449 }
1450
1451 if (h->GetEntries() == 0.0)
1452 {
1453 std::cout << "WARNING! No entries in hcal (" << i << "," << j << "). Skipping this tower" << std::endl;
1454 continue;
1455 }
1456
1457 h_cln = (TH1 *) h->Clone();
1458
1459 h_cln_nm = std::format("h_cln_eta_{}_phi_{}", i, j);
1460
1461 h_cln->SetName(h_cln_nm.c_str());
1462
1463 if (i == 0)
1464 {
1465 scale = pow(10, power - i);
1466 }
1467 else
1468 {
1469 scale = pow(10, power - i - cntr);
1470 }
1471
1472 h_cln->Scale(h_es_cln->Integral(h_es_cln->FindBin(fitMin), h_es_cln->FindBin(fitMax)) / h_cln->Integral(h_cln->FindBin(fitMin), h_cln->FindBin(fitMax)));
1473 h_cln->Scale(scale);
1474 h_cln->GetXaxis()->SetRangeUser(0., xaxisRange);
1475 h_cln->GetYaxis()->SetRangeUser(10, 1e35);
1476 h_cln->Draw("same hist");
1477
1478 if ((i == 0 || i == 12) && j == 0)
1479 {
1480 t->AddEntry(h_cln, "Tower", "l");
1481 }
1482
1483 h_cln = nullptr;
1484
1485 }
1486
1487
1488
1489 h_es_cln->Scale(scale);
1490 h_es_cln->GetXaxis()->SetRangeUser(0., xaxisRange);
1491 h_es_cln->GetYaxis()->SetRangeUser(10, 1e35);
1492 h_es_cln->SetLineColor(2);
1493 h_es_cln->SetLineWidth(3);
1494 h_es_cln->Draw("same hist");
1495
1496 if (i == 0 || i == 12)
1497 {
1498 t->AddEntry(h_es_cln, "#Sigma towers", "l");
1499 }
1500
1501 h_es_cln = nullptr;
1502
1503 if (i > 0)
1504 {
1505 cntr++;
1506 }
1507
1508 }
1509
1510 t->Draw("same");
1511 t = nullptr;
1512 delete t;
1513
1514 }
1515
1516 c->Write();
1517 c = nullptr;
1518 delete c;
1519
1520 }
1521
1522 if (calotype == LiteCaloEval::CEMC)
1523 {
1524 std::cout << "Drawing emcal tower spectra" << std::endl;
1525
1526 int starteta = 0;
1527 int maxeta = 12;
1528
1529 for (int k = 1; k < 9; k++)
1530 {
1531 int power = 36;
1532 int cntr = 3;
1533
1534 if (k > 1)
1535 {
1536 starteta = maxeta;
1537 maxeta = maxeta + 12;
1538 }
1539
1540 TCanvas *c = new TCanvas();
1541 c->SetName(std::format("emcal_eta{}_{}", starteta, (maxeta - 1)).c_str());
1542 gPad->SetLogy(1);
1543
1544 TLegend *t = new TLegend(0.5, 0.8, 0.75, 0.9);
1545 t->SetFillStyle(0);
1546 t->SetBorderSize(0);
1547 t->AddEntry("", std::format("ieta {} - {}", starteta, (maxeta - 1)).c_str());
1548
1549
1550 for (int i = starteta; i < maxeta; i++)
1551 {
1552 etaSliceName = "eta_" + std::to_string(i);
1553
1554 fout->GetObject(etaSliceName.c_str(), h_etaSlice);
1555
1556 if (!h_etaSlice)
1557 {
1558 std::cout << "ERROR! Could not get emcal eta slice " << i << ". Exiting analysis." << std::endl;
1559 exit(-1);
1560 }
1561
1562 if (h_etaSlice->GetEntries() == 0.0)
1563 {
1564 std::cout << "WARNING! EMCal eta slice " << i << " has no entries!" << std::endl;
1565 continue;
1566 }
1567
1568 h_es_cln = (TH1 *) h_etaSlice->Clone();
1569
1570 h_es_nm = "cln_es_" + std::to_string(i);
1571
1572 h_es_cln->SetName(h_es_nm.c_str());
1573
1574
1575 for (int j = 0; j < 256; j++)
1576 {
1577 histName = "emc_ieta" + std::to_string(i) + "_phi" + std::to_string(j);
1578
1579 fout->GetObject(histName.c_str(), h);
1580
1581 if (!h)
1582 {
1583 std::cout << "ERROR! Could not find " << histName << ". Exiting macro." << std::endl;
1584 exit(-1);
1585 }
1586 if (h->GetEntries() == 0.0)
1587 {
1588 std::cout << "WARNING! No entries in emcal (" << i << "," << j << "). Skipping this tower" << std::endl;
1589 continue;
1590 }
1591
1592 h_cln = (TH1 *) h->Clone();
1593
1594 h_cln_nm = std::format("h_cln_eta_{}_phi_{}", i, j);
1595
1596 h_cln->SetName(h_cln_nm.c_str());
1597
1598 if (i == 0 || i % 12 == 0)
1599 {
1600 scale = pow(10, power);
1601 }
1602 else
1603 {
1604 scale = pow(10, power - cntr);
1605 }
1606
1607 h_cln->Scale(h_es_cln->Integral(h_es_cln->FindBin(fitMin), h_es_cln->FindBin(fitMax)) / h_cln->Integral(h_cln->FindBin(fitMin), h_cln->FindBin(fitMax)));
1608 h_cln->Scale(scale);
1609 h_cln->GetXaxis()->SetRangeUser(0., xaxisRange);
1610 h_cln->GetYaxis()->SetRangeUser(0.00001, 1e40);
1611 h_cln->Draw("same hist");
1612
1613 if (((i + 1) % 12 == 0) && j == 0)
1614 {
1615 t->AddEntry(h_cln, "Tower", "l");
1616 }
1617
1618 h_cln = nullptr;
1619
1620 }
1621
1622 h_es_cln->Scale(scale);
1623 h_es_cln->GetXaxis()->SetRangeUser(0., xaxisRange);
1624 h_es_cln->GetYaxis()->SetRangeUser(0.00001, 1e40);
1625 h_es_cln->SetLineColor(2);
1626 h_es_cln->SetLineWidth(3);
1627 h_es_cln->Draw("same hist");
1628
1629 if ((i + 1) % 12 == 0)
1630 {
1631 t->AddEntry(h_es_cln, "#Sigma towers", "l");
1632 }
1633
1634 h_es_cln = nullptr;
1635
1636 if (!(i % 12 == 0))
1637 {
1638 cntr += 3;
1639 }
1640
1641 }
1642
1643 t->Draw("same");
1644
1645 c->Write();
1646
1647 t = nullptr;
1648 delete t;
1649
1650 c = nullptr;
1651 delete c;
1652
1653 std::cout << "Drew canvas " << k << std::endl;
1654 }
1655
1656 }
1657
1658 fout->Close();
1659 fout = nullptr;
1660 delete fout;
1661
1662 std::cout << "Drawing histos is complete." << std::endl;
1663
1664 }
1665
1666 void LiteCaloEval::draw_spectra(const char *infile, const char *outfile)
1667 {
1668 if (calotype == LiteCaloEval::NONE)
1669 {
1670 std::cout << "Did not enter correct calotype. Exiting macro..." << std::endl;
1671 exit(-1);
1672 }
1673
1674 TFile *fin = new TFile(infile, "READ");
1675 TFile *fout = new TFile(outfile, "UPDATE");
1676
1677 TH1F *h = nullptr;
1678 std::string histName;
1679
1680 TH1F *h_etaSlice = nullptr;
1681 std::string etaSliceName;
1682
1683 TH1F *h_cln = nullptr;
1684 std::string h_cln_nm;
1685
1686 TH1F *h_es_cln = nullptr;
1687 std::string h_es_nm;
1688
1689 float scale = 1.0;
1690
1691 float fitMin = getFitMin();
1692 float fitMax = getFitMax();
1693
1694 double xaxisRange = 2.0;
1695
1696 if (calotype == LiteCaloEval::HCALIN)
1697 {
1698 xaxisRange = 3.0;
1699 }
1700
1701 if (calotype == LiteCaloEval::HCALIN || calotype == LiteCaloEval::HCALOUT)
1702 {
1703 std::cout << "Drawing hcal tower spectra" << std::endl;
1704
1705 int starteta = 0;
1706 int maxeta = 12;
1707 int power = 23;
1708
1709 TCanvas *c = new TCanvas();
1710 c->Divide(2);
1711 c->SetName("hcal_spectra");
1712
1713 int cntr = 1;
1714
1715
1716 for (int k = 1; k < 3; k++)
1717 {
1718 c->cd(k);
1719 gPad->SetLogy(1);
1720
1721 if (k == 2)
1722 {
1723 power = 47;
1724 starteta = 12;
1725 maxeta = 24;
1726 }
1727
1728 TLegend *t = new TLegend(0.5, 0.8, 0.75, 0.9);
1729 t->SetFillStyle(0);
1730 t->SetBorderSize(0);
1731 t->AddEntry("", std::format("ieta {} - {}", starteta, (maxeta - 1)).c_str(), "");
1732
1733
1734 for (int i = starteta; i < maxeta; i++)
1735 {
1736
1737 if (calotype == LiteCaloEval::HCALOUT)
1738 {
1739 etaSliceName = "hcalout_eta_" + std::to_string(i);
1740 }
1741 else if (calotype == LiteCaloEval::HCALIN)
1742 {
1743 etaSliceName = "hcalin_eta_" + std::to_string(i);
1744 }
1745
1746 h_etaSlice = (TH1F *) fin->Get(etaSliceName.c_str());
1747
1748 if (!h_etaSlice)
1749 {
1750 std::cout << "ERROR! Could not get hcal eta slice histogram " << i << " ." << std::endl;
1751 gSystem->Exit(1);
1752 }
1753
1754 if (h_etaSlice->GetEntries() == 0.0)
1755 {
1756 std::cout << "WARNING! hcal eta slice " << i << " has no entries!" << std::endl;
1757 continue;
1758 }
1759
1760
1761 double esBW = h_etaSlice->GetBinWidth(1);
1762
1763 double bW = get_spectra_binWidth();
1764
1765 int rbn = std::ceil(bW / esBW);
1766
1767 if (esBW < bW)
1768 {
1769 h_etaSlice->Rebin(rbn);
1770 }
1771
1772 h_es_cln = (TH1F *) h_etaSlice->Clone();
1773
1774 h_es_nm = "cln_es_" + std::to_string(i);
1775
1776 h_es_cln->SetName(h_es_nm.c_str());
1777
1778
1779 for (int j = 0; j < 64; j++)
1780 {
1781 if (calotype == LiteCaloEval::HCALOUT)
1782 {
1783 histName = "hcal_out_eta_" + std::to_string(i) + "_phi_" + std::to_string(j);
1784 }
1785 else if (calotype == LiteCaloEval::HCALIN)
1786 {
1787 histName = "hcal_in_eta_" + std::to_string(i) + "_phi_" + std::to_string(j);
1788 }
1789
1790 h = (TH1F *) fin->Get(histName.c_str());
1791
1792 if (!h)
1793 {
1794 std::cout << "ERROR! Could not find tower " << histName << "." << std::endl;
1795 gSystem->Exit(1);
1796 }
1797
1798 if (h->GetEntries() == 0.0)
1799 {
1800 std::cout << "WARNING! No entries in hcal (" << i << "," << j << "). Skipping this tower" << std::endl;
1801 continue;
1802 }
1803
1804
1805 double tBW = h->GetBinWidth(1);
1806 double W = get_spectra_binWidth();
1807 int rbn2 = std::ceil(W / tBW);
1808
1809 if (tBW < W)
1810 {
1811 h->Rebin(rbn2);
1812 }
1813
1814 h_cln = (TH1F *) h->Clone();
1815
1816 h_cln_nm = std::format("h_cln_eta_{}_phi_{}", i, j);
1817
1818 h_cln->SetName(h_cln_nm.c_str());
1819
1820 if (i == 0)
1821 {
1822 scale = pow(10, power - i);
1823 }
1824 else
1825 {
1826 scale = pow(10, power - i - cntr);
1827 }
1828
1829 h_cln->Scale(h_es_cln->Integral(h_es_cln->FindBin(fitMin), h_es_cln->FindBin(fitMax)) / h_cln->Integral(h_cln->FindBin(fitMin), h_cln->FindBin(fitMax)));
1830 h_cln->Scale(scale);
1831 h_cln->GetXaxis()->SetRangeUser(0., xaxisRange);
1832 h_cln->GetYaxis()->SetRangeUser(10, 1e35);
1833 h_cln->Draw("same hist");
1834
1835 if ((i == 0 || i == 12) && j == 0)
1836 {
1837 t->AddEntry(h_cln, "Tower", "l");
1838 }
1839
1840 h_cln = nullptr;
1841
1842 }
1843
1844
1845
1846 h_es_cln->Scale(scale);
1847 h_es_cln->GetXaxis()->SetRangeUser(0., xaxisRange);
1848 h_es_cln->GetYaxis()->SetRangeUser(10, 1e35);
1849 h_es_cln->SetLineColor(2);
1850 h_es_cln->SetLineWidth(3);
1851 h_es_cln->Draw("same hist");
1852
1853 if (i == 0 || i == 12)
1854 {
1855 t->AddEntry(h_es_cln, "#Sigma towers", "l");
1856 }
1857
1858 h_es_cln = nullptr;
1859
1860 if (i > 0)
1861 {
1862 cntr++;
1863 }
1864
1865 }
1866
1867 t->Draw("same");
1868 t = nullptr;
1869 delete t;
1870
1871 }
1872
1873 fout->cd();
1874 c->Write();
1875 c = nullptr;
1876 delete c;
1877
1878 }
1879
1880 if (calotype == LiteCaloEval::CEMC)
1881 {
1882 std::cout << "Drawing emcal tower spectra" << std::endl;
1883
1884 int starteta = 0;
1885 int maxeta = 12;
1886
1887 for (int k = 1; k < 9; k++)
1888 {
1889 int power = 36;
1890 int cntr = 3;
1891
1892 if (k > 1)
1893 {
1894 starteta = maxeta;
1895 maxeta = maxeta + 12;
1896 }
1897
1898 TCanvas *c = new TCanvas();
1899 c->SetName(std::format("emcal_eta{}_{}", starteta, (maxeta - 1)).c_str());
1900 gPad->SetLogy(1);
1901
1902 TLegend *t = new TLegend(0.5, 0.8, 0.75, 0.9);
1903 t->AddEntry("", std::format("ieta {} - {}", starteta, (maxeta - 1)).c_str());
1904
1905
1906 for (int i = starteta; i < maxeta; i++)
1907 {
1908 etaSliceName = "eta_" + std::to_string(i);
1909
1910 h_etaSlice = (TH1F *) fin->Get(etaSliceName.c_str());
1911
1912 if (!h_etaSlice)
1913 {
1914 std::cout << "ERROR! Could not get emcal eta slice " << i << ". Exiting analysis." << std::endl;
1915 exit(-1);
1916 }
1917
1918 if (h_etaSlice->GetEntries() == 0.0)
1919 {
1920 std::cout << "WARNING! EMCal eta slice " << i << " has no entries!" << std::endl;
1921 continue;
1922 }
1923
1924
1925 double esBW = h_etaSlice->GetBinWidth(1);
1926 double bW = get_spectra_binWidth();
1927 int rbn = std::ceil(bW / esBW);
1928
1929 if (esBW < bW)
1930 {
1931 h_etaSlice->Rebin(rbn);
1932 }
1933
1934 h_es_cln = (TH1F *) h_etaSlice->Clone();
1935
1936 h_es_nm = "cln_es_" + std::to_string(i);
1937
1938 h_es_cln->SetName(h_es_nm.c_str());
1939
1940
1941 for (int j = 0; j < 256; j++)
1942 {
1943 histName = "emc_ieta" + std::to_string(i) + "_phi" + std::to_string(j);
1944
1945 h = (TH1F *) fin->Get(histName.c_str());
1946
1947 if (!h)
1948 {
1949 std::cout << "ERROR! Could not find " << histName << ". Exiting macro." << std::endl;
1950 exit(-1);
1951 }
1952 if (h->GetEntries() == 0.0)
1953 {
1954 std::cout << "WARNING! No entries in emcal (" << i << "," << j << "). Skipping this tower" << std::endl;
1955 continue;
1956 }
1957
1958
1959 double tBW = h->GetBinWidth(1);
1960 double W = get_spectra_binWidth();
1961 int rbn2 = std::ceil(W / tBW);
1962
1963 if (tBW < W)
1964 {
1965 h->Rebin(rbn2);
1966 }
1967
1968 h_cln = (TH1F *) h->Clone();
1969
1970 h_cln_nm = std::format("h_cln_eta_{}_phi_{}", i, j);
1971
1972 h_cln->SetName(h_cln_nm.c_str());
1973
1974 if (i == 0 || i % 12 == 0)
1975 {
1976 scale = pow(10, power);
1977 }
1978 else
1979 {
1980 scale = pow(10, power - cntr);
1981 }
1982
1983 h_cln->Scale(h_es_cln->Integral(h_es_cln->FindBin(fitMin), h_es_cln->FindBin(fitMax)) / h_cln->Integral(h_cln->FindBin(fitMin), h_cln->FindBin(fitMax)));
1984 h_cln->Scale(scale);
1985 h_cln->GetXaxis()->SetRangeUser(0., xaxisRange);
1986 h_cln->GetYaxis()->SetRangeUser(0.00001, 1e40);
1987 h_cln->Draw("same hist");
1988
1989 if (((i + 1) % 12 == 0) && j == 0)
1990 {
1991 t->AddEntry(h_cln, "Tower", "l");
1992 }
1993
1994 h_cln = nullptr;
1995
1996 }
1997
1998 h_es_cln->Scale(scale);
1999 h_es_cln->GetXaxis()->SetRangeUser(0., xaxisRange);
2000 h_es_cln->GetYaxis()->SetRangeUser(0.00001, 1e40);
2001 h_es_cln->SetLineColor(2);
2002 h_es_cln->SetLineWidth(3);
2003 h_es_cln->Draw("same hist");
2004
2005 if ((i + 1) % 12 == 0)
2006 {
2007 t->AddEntry(h_es_cln, "#Sigma towers", "l");
2008 }
2009
2010 h_es_cln = nullptr;
2011
2012 if (!(i % 12 == 0))
2013 {
2014 cntr += 3;
2015 }
2016
2017 }
2018
2019 t->Draw("same");
2020
2021 fout->cd();
2022
2023 c->Write();
2024
2025 t = nullptr;
2026 delete t;
2027
2028 c = nullptr;
2029 delete c;
2030
2031 std::cout << "Drew canvas " << k << std::endl;
2032 }
2033
2034 }
2035
2036 fout->Close();
2037 fout = nullptr;
2038 delete fout;
2039
2040 std::cout << "Drawing histos is complete." << std::endl;
2041 }
2042
2043 void LiteCaloEval::fit_info(const char *outfile, const int runNum)
2044 {
2045 TFile *fout = new TFile(outfile, "UPDATE");
2046
2047 int eta;
2048 int phi;
2049 double towers;
2050
2051
2052 if (calotype == LiteCaloEval::HCALIN || calotype == LiteCaloEval::HCALOUT)
2053 {
2054 eta = 24;
2055 phi = 64;
2056 towers = 1536.0;
2057 }
2058 else if (calotype == LiteCaloEval::CEMC)
2059 {
2060 eta = 96;
2061 phi = 256;
2062 towers = 24576.0;
2063 }
2064 else
2065 {
2066 std::cout << "calotype not set. Exiting." << std::endl;
2067 exit(-1);
2068 }
2069
2070 TH2 *errMap = new TH2F("errMap", "", eta, 0, eta, phi, 0, phi);
2071 errMap->GetXaxis()->SetTitle("#eta Bin");
2072 errMap->GetYaxis()->SetTitle("#phi Bin");
2073
2074
2075 TH1 *chi2 = new TH1F("chi2", "", 500, 0, 50);
2076 chi2->GetXaxis()->SetTitle("#chi^{2} / NDF");
2077 chi2->GetYaxis()->SetTitle("Counts");
2078
2079
2080 TH2 *chi2Map = new TH2F("chi2Map", "", eta, 0, eta, phi, 0, phi);
2081 chi2Map->GetXaxis()->SetTitle("#eta Bin");
2082 chi2Map->GetYaxis()->SetTitle("#phi Bin");
2083
2084
2085 TH2 *fitFail = new TH2F("fitFail", "", eta, 0, eta, phi, 0, phi);
2086 fitFail->GetXaxis()->SetTitle("#eta bin");
2087 fitFail->GetYaxis()->SetTitle("#phi bin");
2088
2089 std::string histname;
2090 TH1 *htmp = nullptr;
2091 TF1 *fn = nullptr;
2092
2093
2094 TH2 *cp = nullptr;
2095 fout->GetObject("corrPat", cp);
2096
2097 if (!cp)
2098 {
2099 std::cout << "Error! Did not get corrPat histogram. Exiting analysis." << std::endl;
2100 exit(-1);
2101 }
2102
2103 TGraphErrors *tmp_avgTSC = new TGraphErrors();
2104 tmp_avgTSC->SetName("g_avgTSC");
2105 tmp_avgTSC->GetXaxis()->SetTitle("run number");
2106 tmp_avgTSC->GetYaxis()->SetTitle("Mean Towerslope Correction");
2107 tmp_avgTSC->SetMarkerStyle(8);
2108 tmp_avgTSC->SetMarkerSize(1);
2109
2110 double sum4avg = 0.0;
2111 double tscAvg = 0.0;
2112 double sum4SE = 0.0;
2113 double SE = 0.0;
2114
2115
2116 for (int i = 0; i < eta; i++)
2117 {
2118
2119 for (int j = 0; j < phi; j++)
2120 {
2121
2122 if (calotype == LiteCaloEval::HCALOUT)
2123 {
2124 histname = std::format("hcal_out_eta_{}_phi_{}", i, j);
2125 }
2126
2127
2128 if (calotype == LiteCaloEval::HCALIN)
2129 {
2130 histname = std::format("hcal_in_eta_{}_phi_{}", i, j);
2131 }
2132
2133
2134 if (calotype == LiteCaloEval::CEMC)
2135 {
2136 histname = std::format("emc_ieta{}_phi{}", i, j);
2137 }
2138
2139 sum4avg += (cp->GetBinContent(i + 1, j + 1));
2140
2141 fout->GetObject(histname.c_str(), htmp);
2142
2143 fn = htmp->GetFunction("myexpo");
2144
2145 errMap->SetBinContent(i + 1, j + 1, fn->GetParError(1));
2146
2147 chi2->Fill(fn->GetChisquare() / fn->GetNDF());
2148
2149 chi2Map->SetBinContent(i + 1, j + 1, fn->GetChisquare() / fn->GetNDF());
2150
2151 if (fn->GetChisquare() / fn->GetNDF() > 5)
2152 {
2153 fitFail->Fill(i, j);
2154
2155 }
2156
2157 }
2158
2159 }
2160
2161 tscAvg = sum4avg / towers;
2162
2163
2164 for (int i = 0; i < eta; i++)
2165 {
2166 for (int j = 0; j < phi; j++)
2167 {
2168 sum4SE += pow(cp->GetBinContent(i + 1, j + 1) - tscAvg, 2);
2169 }
2170 }
2171
2172 SE = sqrt(sum4SE / towers) / sqrt(towers);
2173
2174 tmp_avgTSC->SetPoint(0, runNum, tscAvg);
2175 tmp_avgTSC->SetPointError(0, 0.0, SE);
2176
2177 errMap->Write();
2178 chi2->Write();
2179 chi2Map->Write();
2180 fitFail->Write();
2181 tmp_avgTSC->Write();
2182
2183 errMap = nullptr;
2184 delete errMap;
2185
2186 chi2 = nullptr;
2187 delete chi2;
2188
2189 chi2Map = nullptr;
2190 delete chi2Map;
2191
2192 fitFail = nullptr;
2193 delete fitFail;
2194
2195 fout->Close();
2196 fout = nullptr;
2197 delete fout;
2198
2199 std::cout << "Finished fit info" << std::endl;
2200 }
2201
2202
2203 void LiteCaloEval::fit_info(const char *infile, const char *outfile, const int runNum)
2204 {
2205 TFile *fin = new TFile(infile, "READ");
2206 TFile *fout = new TFile(outfile, "UPDATE");
2207
2208 int eta;
2209 int phi;
2210 double towers;
2211
2212 if (calotype == LiteCaloEval::HCALIN || calotype == LiteCaloEval::HCALOUT)
2213 {
2214 eta = 24;
2215 phi = 64;
2216 towers = 1536.0;
2217 }
2218 else if (calotype == LiteCaloEval::CEMC)
2219 {
2220 eta = 96;
2221 phi = 256;
2222 towers = 24576.0;
2223 }
2224 else
2225 {
2226 std::cout << "calotype not set. Exiting." << std::endl;
2227 exit(-1);
2228 }
2229
2230 TH2 *errMap = new TH2F("errMap", "", eta, 0, eta, phi, 0, phi);
2231 errMap->GetXaxis()->SetTitle("#eta Bin");
2232 errMap->GetYaxis()->SetTitle("#phi Bin");
2233
2234
2235 TH1 *chi2 = new TH1F("chi2", "", 500, 0, 50);
2236 chi2->GetXaxis()->SetTitle("#chi^{2} / NDF");
2237 chi2->GetYaxis()->SetTitle("Counts");
2238
2239
2240 TH2 *chi2Map = new TH2F("chi2Map", "", eta, 0, eta, phi, 0, phi);
2241 chi2Map->GetXaxis()->SetTitle("#eta Bin");
2242 chi2Map->GetYaxis()->SetTitle("#phi Bin");
2243
2244
2245 TH2 *fitFail = new TH2F("fitFail", "", eta, 0, eta, phi, 0, phi);
2246 fitFail->GetXaxis()->SetTitle("#eta bin");
2247 fitFail->GetYaxis()->SetTitle("#phi bin");
2248
2249 std::string histname;
2250 TH1 *htmp{nullptr};
2251 TF1 *fn{nullptr};
2252
2253
2254 TH2 *cp{nullptr};
2255 fin->GetObject("corrPat", cp);
2256
2257 if (!cp)
2258 {
2259 std::cout << "Error! Did not get corrPat histogram. Exiting analysis." << std::endl;
2260 exit(-1);
2261 }
2262
2263 TGraphErrors *tmp_avgTSC = new TGraphErrors();
2264 tmp_avgTSC->SetName("g_avgTSC");
2265 tmp_avgTSC->GetXaxis()->SetTitle("run number");
2266 tmp_avgTSC->GetYaxis()->SetTitle("Mean Towerslope Correction");
2267 tmp_avgTSC->SetMarkerStyle(8);
2268 tmp_avgTSC->SetMarkerSize(1);
2269
2270 double sum4avg = 0.0;
2271 double tscAvg = 0.0;
2272 double sum4SE = 0.0;
2273 double SE = 0.0;
2274
2275
2276 for (int i = 0; i < eta; i++)
2277 {
2278
2279 for (int j = 0; j < phi; j++)
2280 {
2281
2282 if (calotype == LiteCaloEval::HCALOUT)
2283 {
2284 histname = std::format("hcal_out_eta_{}_phi_{}", i, j);
2285 }
2286
2287
2288 if (calotype == LiteCaloEval::HCALIN)
2289 {
2290 histname = std::format("hcal_in_eta_{}_phi_{}", i, j);
2291 }
2292
2293
2294 if (calotype == LiteCaloEval::CEMC)
2295 {
2296 histname = std::format("emc_ieta{}_phi{}", i, j);
2297 }
2298
2299 sum4avg += (cp->GetBinContent(i + 1, j + 1));
2300
2301 fin->GetObject(histname.c_str(), htmp);
2302
2303 fn = htmp->GetFunction("myexpo");
2304
2305 errMap->SetBinContent(i + 1, j + 1, fn->GetParError(1));
2306
2307 chi2->Fill(fn->GetChisquare() / fn->GetNDF());
2308
2309 chi2Map->SetBinContent(i + 1, j + 1, fn->GetChisquare() / fn->GetNDF());
2310
2311 if (fn->GetChisquare() / fn->GetNDF() > 5)
2312 {
2313 fitFail->Fill(i, j);
2314 }
2315
2316 }
2317
2318 }
2319
2320 tscAvg = sum4avg / towers;
2321
2322
2323 for (int i = 0; i < eta; i++)
2324 {
2325 for (int j = 0; j < phi; j++)
2326 {
2327 sum4SE += pow(cp->GetBinContent(i + 1, j + 1) - tscAvg, 2);
2328 }
2329 }
2330
2331 SE = sqrt(sum4SE / towers) / sqrt(towers);
2332
2333 tmp_avgTSC->SetPoint(0, runNum, tscAvg);
2334 tmp_avgTSC->SetPointError(0, 0.0, SE);
2335
2336 errMap->Write();
2337 chi2->Write();
2338 chi2Map->Write();
2339 fitFail->Write();
2340 tmp_avgTSC->Write();
2341
2342 errMap = nullptr;
2343 delete errMap;
2344
2345 chi2 = nullptr;
2346 delete chi2;
2347
2348 chi2Map = nullptr;
2349 delete chi2Map;
2350
2351 fitFail = nullptr;
2352 delete fitFail;
2353
2354 fin->Close();
2355 delete fin;
2356 fin = nullptr;
2357
2358 fout->Close();
2359 fout = nullptr;
2360 delete fout;
2361
2362 std::cout << "Finished fit info" << std::endl;
2363 }