File indexing completed on 2025-08-05 08:11:14
0001 #include "TriggerPerformance.h"
0002 #include "FindDijets.h"
0003
0004 #include <TH1.h>
0005 #include <TH2.h>
0006 #include <TH2D.h>
0007 #include <TFile.h>
0008 #include <TF1.h>
0009 #include <TH1D.h>
0010 #include <g4main/PHG4TruthInfoContainer.h>
0011 #include <g4main/PHG4Particle.h>
0012 #include <TMath.h>
0013 #include <calotrigger/TriggerRunInfov1.h>
0014
0015 #include <calobase/RawTowerGeom.h>
0016 #include <jetbackground/TowerBackground.h>
0017 #include <jetbase/Jet.h>
0018 #include <jetbase/JetContainer.h>
0019 #include <jetbase/JetContainerv1.h>
0020 #include <jetbase/Jetv2.h>
0021
0022 #include <calobase/TowerInfoContainer.h>
0023 #include <calobase/TowerInfoContainerv1.h>
0024 #include <calobase/TowerInfoContainerv2.h>
0025 #include <calobase/TowerInfoContainerv3.h>
0026 #include <calobase/TowerInfo.h>
0027 #include <calobase/TowerInfov1.h>
0028 #include <calobase/TowerInfov2.h>
0029 #include <calobase/TowerInfov3.h>
0030
0031 #include <calobase/RawCluster.h>
0032 #include <calobase/RawClusterContainer.h>
0033 #include <calobase/RawClusterUtility.h>
0034 #include <calobase/RawTowerGeomContainer.h>
0035 #include <calobase/RawTower.h>
0036 #include <calobase/RawTowerContainer.h>
0037 #include <fun4all/Fun4AllHistoManager.h>
0038 #include <HepMC/SimpleVector.h>
0039
0040 #include <globalvertex/GlobalVertex.h>
0041 #include <globalvertex/GlobalVertexMap.h>
0042 #include "DijetEventDisplay.h"
0043 #include <vector>
0044
0045 #include <fun4all/Fun4AllReturnCodes.h>
0046 #include <phool/PHCompositeNode.h>
0047 #include <phool/PHIODataNode.h>
0048 #include <phool/PHNode.h>
0049 #include <phool/PHNodeIterator.h>
0050 #include <phool/PHObject.h>
0051 #include <phool/getClass.h>
0052
0053
0054 #include <iostream>
0055
0056 #include <map>
0057
0058
0059 TriggerPerformance::TriggerPerformance(const std::string &name, const std::string &outfilename):
0060 SubsysReco(name)
0061
0062 {
0063 _foutname = outfilename;
0064 m_calo_nodename = "TOWERINFO_CALIB";
0065 m_photon_triggernames = {"Photon 2 GeV+ MBD NS >= 1",
0066 "Photon 3 GeV + MBD NS >= 1",
0067 "Photon 4 GeV + MBD NS >= 1",
0068 "Photon 5 GeV + MBD NS >= 1"};
0069 m_photon_prescales = {-1,
0070 -1,
0071 -1,
0072 -1};
0073
0074 m_jet_triggernames = {"Jet 6 GeV + MBD NS >= 1",
0075 "Jet 8 GeV + MBD NS >= 1",
0076 "Jet 10 GeV + MBD NS >= 1",
0077 "Jet 12 GeV + MBD NS >= 1"};
0078 m_jet_prescales = {-1,
0079 -1,
0080 -1,
0081 -1};
0082
0083
0084 m_jet_emu_triggernames = {"Jet 6 GeV",
0085 "Jet 8 GeV",
0086 "Jet 10 GeV",
0087 "Jet 12 GeV"};
0088 m_photon_emu_triggernames = {"Photon 2 GeV",
0089 "Photon 3 GeV",
0090 "Photon 4 GeV",
0091 "Photon 5 GeV"};
0092
0093 m_mbd_prescale = -1;
0094 }
0095
0096
0097 TriggerPerformance::~TriggerPerformance()
0098 {
0099
0100 }
0101
0102
0103 int TriggerPerformance::Init(PHCompositeNode *topNode)
0104 {
0105 _i_event = 0;
0106
0107
0108 triggeranalyzer = new TriggerAnalyzer();
0109
0110 triggeranalyzer->UseEmulator(m_useEmulator);
0111
0112 hm = new Fun4AllHistoManager("TRIGGER_HISTOGRAMS");
0113
0114 TH1D *h;
0115 TH2D *h2;
0116
0117
0118 h = new TH1D("h_jet_et_mb", ";E_{T}; count", 200, 0, 100);
0119 hm->registerHisto(h);
0120 h = new TH1D("h_jet_pt_mb", ";p_{T}; count", 200, 0, 100);
0121 hm->registerHisto(h);
0122
0123 h = new TH1D("h_jet_et_scaled_mb", ";E_{T}; count", 200, 0, 100);
0124 hm->registerHisto(h);
0125 h = new TH1D("h_jet_pt_scaled_mb", ";p_{T}; count", 200, 0, 100);
0126 hm->registerHisto(h);
0127
0128 h2 = new TH2D("h_2d_jet_dist_mb",";#eta;#phi", 24, -1.1, 1.1, 64, -1*TMath::Pi(), TMath::Pi());
0129 hm->registerHisto(h2);
0130
0131 h = new TH1D("h_jet_max_et_mb", ";E_{T}; count", 200, 0, 100);
0132 hm->registerHisto(h);
0133 h = new TH1D("h_jet_max_pt_mb", ";p_{T}; count", 200, 0, 100);
0134 hm->registerHisto(h);
0135
0136 h = new TH1D("h_jet_max_et_scaled_mb", ";E_{T}; count", 200, 0, 100);
0137 hm->registerHisto(h);
0138 h = new TH1D("h_jet_max_pt_scaled_mb", ";p_{T}; count", 200, 0, 100);
0139 hm->registerHisto(h);
0140
0141 h2 = new TH2D("h_2d_jet_max_dist_mb",";#eta;#phi", 24, -1.1, 1.1, 64, -1*TMath::Pi(), TMath::Pi());
0142 hm->registerHisto(h2);
0143
0144 for (int j = 0; j < 4; j++)
0145 {
0146 h2 = new TH2D(Form("h_2d_jet_dist_trig_%d", j),";#eta;#phi", 24, -1.1, 1.1, 64, -1*TMath::Pi(), TMath::Pi());
0147 hm->registerHisto(h2);
0148
0149
0150 h = new TH1D(Form("h_jet_et_trig_%d", j), ";E_{T}; count", 200, 0, 100);
0151 hm->registerHisto(h);
0152 h = new TH1D(Form("h_jet_pt_trig_%d", j), ";p_{T}; count", 200, 0, 100);
0153 hm->registerHisto(h);
0154
0155 h = new TH1D(Form("h_jet_et_scaled_trig_%d", j), ";E_{T}; count", 200, 0, 100);
0156 hm->registerHisto(h);
0157 h = new TH1D(Form("h_jet_pt_scaled_trig_%d", j), ";p_{T}; count", 200, 0, 100);
0158 hm->registerHisto(h);
0159
0160 h = new TH1D(Form("h_jet_et_raw_trig_%d", j), ";E_{T}; count", 200, 0, 100);
0161 hm->registerHisto(h);
0162 h = new TH1D(Form("h_jet_pt_raw_trig_%d", j), ";p_{T}; count", 200, 0, 100);
0163 hm->registerHisto(h);
0164
0165 h2 = new TH2D(Form("h_2d_jet_max_dist_trig_%d", j),";#eta;#phi", 24, -1.1, 1.1, 64, -1*TMath::Pi(), TMath::Pi());
0166 hm->registerHisto(h2);
0167
0168
0169 h = new TH1D(Form("h_jet_max_et_trig_%d", j), ";E_{T}; count", 200, 0, 100);
0170 hm->registerHisto(h);
0171 h = new TH1D(Form("h_jet_max_pt_trig_%d", j), ";p_{T}; count", 200, 0, 100);
0172 hm->registerHisto(h);
0173
0174 h = new TH1D(Form("h_jet_max_et_scaled_trig_%d", j), ";E_{T}; count", 200, 0, 100);
0175 hm->registerHisto(h);
0176 h = new TH1D(Form("h_jet_max_pt_scaled_trig_%d", j), ";p_{T}; count", 200, 0, 100);
0177 hm->registerHisto(h);
0178
0179 h = new TH1D(Form("h_jet_max_et_raw_trig_%d", j), ";E_{T}; count", 200, 0, 100);
0180 hm->registerHisto(h);
0181 h = new TH1D(Form("h_jet_max_pt_raw_trig_%d", j), ";p_{T}; count", 200, 0, 100);
0182 hm->registerHisto(h);
0183
0184 }
0185 h = new TH1D("h_cluster_et_mb", ";E_{T}; count", 200, 0, 50);
0186 hm->registerHisto(h);
0187 h = new TH1D("h_cluster_pt_mb", ";p_{T}; count", 200, 0, 50);
0188 hm->registerHisto(h);
0189
0190 h = new TH1D("h_cluster_et_scaled_mb", ";E_{T}; count", 200, 0, 50);
0191 hm->registerHisto(h);
0192 h = new TH1D("h_cluster_pt_scaled_mb", ";p_{T}; count", 200, 0, 50);
0193 hm->registerHisto(h);
0194
0195 h2 = new TH2D("h_2d_cluster_dist_mb",";#eta;#phi", 24, -1.1, 1.1, 64, -1*TMath::Pi(), TMath::Pi());
0196 hm->registerHisto(h2);
0197
0198 h = new TH1D("h_cluster_max_et_mb", ";E_{T}; count", 200, 0, 50);
0199 hm->registerHisto(h);
0200 h = new TH1D("h_cluster_max_pt_mb", ";p_{T}; count", 200, 0, 50);
0201 hm->registerHisto(h);
0202
0203 h = new TH1D("h_cluster_max_et_scaled_mb", ";E_{T}; count", 200, 0, 50);
0204 hm->registerHisto(h);
0205 h = new TH1D("h_cluster_max_pt_scaled_mb", ";p_{T}; count", 200, 0, 50);
0206 hm->registerHisto(h);
0207
0208 h2 = new TH2D("h_2d_cluster_max_dist_mb",";#eta;#phi", 24, -1.1, 1.1, 64, -1*TMath::Pi(), TMath::Pi());
0209 hm->registerHisto(h2);
0210
0211 for (int j = 0; j < 4; j++)
0212 {
0213 h2 = new TH2D(Form("h_2d_cluster_dist_trig_%d", j),";#eta;#phi", 24, -1.1, 1.1, 64, -1*TMath::Pi(), TMath::Pi());
0214 hm->registerHisto(h2);
0215
0216
0217 h = new TH1D(Form("h_cluster_et_trig_%d", j), ";E_{T}; count", 200, 0, 50);
0218 hm->registerHisto(h);
0219 h = new TH1D(Form("h_cluster_pt_trig_%d", j), ";p_{T}; count", 200, 0, 50);
0220 hm->registerHisto(h);
0221
0222 h = new TH1D(Form("h_cluster_et_scaled_trig_%d", j), ";E_{T}; count", 200, 0, 50);
0223 hm->registerHisto(h);
0224 h = new TH1D(Form("h_cluster_pt_scaled_trig_%d", j), ";p_{T}; count", 200, 0, 50);
0225 hm->registerHisto(h);
0226
0227 h = new TH1D(Form("h_cluster_et_raw_trig_%d", j), ";E_{T}; count", 200, 0, 50);
0228 hm->registerHisto(h);
0229 h = new TH1D(Form("h_cluster_pt_raw_trig_%d", j), ";p_{T}; count", 200, 0, 50);
0230 hm->registerHisto(h);
0231
0232 h2 = new TH2D(Form("h_2d_cluster_max_dist_trig_%d", j),";#eta;#phi", 24, -1.1, 1.1, 64, -1*TMath::Pi(), TMath::Pi());
0233 hm->registerHisto(h2);
0234
0235
0236 h = new TH1D(Form("h_cluster_max_et_trig_%d", j), ";E_{T}; count", 200, 0, 50);
0237 hm->registerHisto(h);
0238 h = new TH1D(Form("h_cluster_max_pt_trig_%d", j), ";p_{T}; count", 200, 0, 50);
0239 hm->registerHisto(h);
0240
0241 h = new TH1D(Form("h_cluster_max_et_scaled_trig_%d", j), ";E_{T}; count", 200, 0, 50);
0242 hm->registerHisto(h);
0243 h = new TH1D(Form("h_cluster_max_pt_scaled_trig_%d", j), ";p_{T}; count", 200, 0, 50);
0244 hm->registerHisto(h);
0245
0246 h = new TH1D(Form("h_cluster_max_et_raw_trig_%d", j), ";E_{T}; count", 200, 0, 50);
0247 hm->registerHisto(h);
0248 h = new TH1D(Form("h_cluster_max_pt_raw_trig_%d", j), ";p_{T}; count", 200, 0, 50);
0249 hm->registerHisto(h);
0250
0251 }
0252
0253 return Fun4AllReturnCodes::EVENT_OK;
0254 }
0255
0256
0257 int TriggerPerformance::InitRun(PHCompositeNode *topNode)
0258 {
0259 if (m_useEmulator)
0260 {
0261 std::cout << "using emulator " <<std::endl;
0262 m_mbd_prescale = 1;
0263 }
0264 else
0265 {
0266 triggeranalyzer->decodeTriggers(topNode);
0267 m_mbd_prescale = triggeranalyzer->getTriggerPrescale("MBD N&S >= 1") + 1;
0268 }
0269 return Fun4AllReturnCodes::EVENT_OK;
0270
0271 }
0272
0273 int TriggerPerformance::process_event(PHCompositeNode *topNode)
0274 {
0275
0276
0277 _i_event++;
0278
0279 triggeranalyzer->decodeTriggers(topNode);
0280 if (Verbosity())
0281 {
0282 triggeranalyzer->Print();
0283 }
0284 std::string recoJetName = "AntiKt_Tower_r04_Sub1";
0285
0286 JetContainer *jets_4 = findNode::getClass<JetContainer>(topNode, recoJetName);
0287
0288 if (!jets_4)
0289 {
0290 std::cout << " Nothing in the jets"<<std::endl;
0291 return Fun4AllReturnCodes::ABORTRUN;
0292 }
0293 GlobalVertexMap* vertexmap = findNode::getClass<GlobalVertexMap>(topNode, "GlobalVertexMap");
0294
0295 float vtx_z = 999;
0296 if (vertexmap && !vertexmap->empty())
0297 {
0298 GlobalVertex* vtx = vertexmap->begin()->second;
0299 if (vtx)
0300 {
0301 vtx_z = vtx->get_z();
0302 }
0303 }
0304 if (vtx_z == 999)
0305 {
0306 return Fun4AllReturnCodes::EVENT_OK;
0307 }
0308
0309 RawClusterContainer *clusters = findNode::getClass<RawClusterContainer>(topNode, "CLUSTERINFO_CEMC");
0310
0311 RawTowerGeomContainer *tower_geomIH = findNode::getClass<RawTowerGeomContainer>(topNode, "TOWERGEOM_HCALIN");
0312 RawTowerGeomContainer *tower_geomOH = findNode::getClass<RawTowerGeomContainer>(topNode, "TOWERGEOM_HCALOUT");
0313
0314 TowerInfoContainer *hcalin_towers = findNode::getClass<TowerInfoContainer>(topNode, m_calo_nodename + "_HCALIN");
0315 TowerInfoContainer *hcalout_towers = findNode::getClass<TowerInfoContainer>(topNode, m_calo_nodename + "_HCALOUT");
0316 TowerInfoContainer *emcalre_towers = findNode::getClass<TowerInfoContainer>(topNode, m_calo_nodename + "_CEMC_RETOWER");
0317
0318 float background_v2 = 0;
0319 float background_Psi2 = 0;
0320 bool has_tower_background = false;
0321 TowerBackground *towBack = findNode::getClass<TowerBackground>(topNode, "TowerInfoBackground_Sub2");
0322 if (towBack)
0323 {
0324 has_tower_background = true;
0325 background_v2 = towBack->get_v2();
0326 background_Psi2 = towBack->get_Psi2();
0327 }
0328
0329 std::vector<jetty> jets{};
0330 for (auto jet : *jets_4)
0331 {
0332 float jet_total_eT = 0;
0333 float eFrac_ihcal = 0;
0334 float eFrac_ohcal = 0;
0335 float eFrac_emcal = 0;
0336
0337 for (auto comp : jet->get_comp_vec())
0338 {
0339
0340 unsigned int channel = comp.second;
0341 TowerInfo *tower;
0342 float tower_eT = 0;
0343 if (comp.first == 26 || comp.first == 30)
0344 {
0345 tower = hcalin_towers->get_tower_at_channel(channel);
0346 if (!tower || !tower_geomIH)
0347 {
0348 continue;
0349 }
0350 if(!tower->get_isGood()) continue;
0351
0352 unsigned int calokey = hcalin_towers->encode_key(channel);
0353
0354 int ieta = hcalin_towers->getTowerEtaBin(calokey);
0355 int iphi = hcalin_towers->getTowerPhiBin(calokey);
0356 const RawTowerDefs::keytype key = RawTowerDefs::encode_towerid(RawTowerDefs::CalorimeterId::HCALIN, ieta, iphi);
0357 float tower_phi = tower_geomIH->get_tower_geometry(key)->get_phi();
0358 float tower_eta = tower_geomIH->get_tower_geometry(key)->get_eta();
0359 tower_eT = tower->get_energy() / std::cosh(tower_eta);
0360
0361 if (comp.first == 30)
0362 {
0363 if (has_tower_background)
0364 {
0365 float UE = towBack->get_UE(1).at(ieta);
0366 float tower_UE = UE * (1 + 2 * background_v2 * std::cos(2 * (tower_phi - background_Psi2)));
0367 tower_eT = (tower->get_energy() - tower_UE) / std::cosh(tower_eta);
0368 }
0369 }
0370
0371 eFrac_ihcal += tower_eT;
0372 jet_total_eT += tower_eT;
0373 }
0374 else if (comp.first == 27 || comp.first == 31)
0375 {
0376 tower = hcalout_towers->get_tower_at_channel(channel);
0377
0378 if (!tower || !tower_geomOH)
0379 {
0380 continue;
0381 }
0382
0383 unsigned int calokey = hcalout_towers->encode_key(channel);
0384 int ieta = hcalout_towers->getTowerEtaBin(calokey);
0385 int iphi = hcalout_towers->getTowerPhiBin(calokey);
0386 const RawTowerDefs::keytype key = RawTowerDefs::encode_towerid(RawTowerDefs::CalorimeterId::HCALOUT, ieta, iphi);
0387 float tower_phi = tower_geomOH->get_tower_geometry(key)->get_phi();
0388 float tower_eta = tower_geomOH->get_tower_geometry(key)->get_eta();
0389 tower_eT = tower->get_energy() / std::cosh(tower_eta);
0390
0391 if (comp.first == 31)
0392 {
0393 if (has_tower_background)
0394 {
0395 float UE = towBack->get_UE(2).at(ieta);
0396 float tower_UE = UE * (1 + 2 * background_v2 * std::cos(2 * (tower_phi - background_Psi2)));
0397 tower_eT = (tower->get_energy() - tower_UE) / std::cosh(tower_eta);
0398 }
0399 }
0400
0401 eFrac_ohcal += tower_eT;
0402 jet_total_eT += tower_eT;
0403 }
0404 else if (comp.first == 28 || comp.first == 29)
0405 {
0406 tower = emcalre_towers->get_tower_at_channel(channel);
0407
0408 if (!tower || !tower_geomIH)
0409 {
0410 continue;
0411 }
0412
0413 unsigned int calokey = emcalre_towers->encode_key(channel);
0414 int ieta = emcalre_towers->getTowerEtaBin(calokey);
0415 int iphi = emcalre_towers->getTowerPhiBin(calokey);
0416 const RawTowerDefs::keytype key = RawTowerDefs::encode_towerid(RawTowerDefs::CalorimeterId::HCALIN, ieta, iphi);
0417 float tower_phi = tower_geomIH->get_tower_geometry(key)->get_phi();
0418 float tower_eta = tower_geomIH->get_tower_geometry(key)->get_eta();
0419 tower_eT = tower->get_energy() / std::cosh(tower_eta);
0420
0421 if (comp.first == 29)
0422 {
0423 if (has_tower_background)
0424 {
0425 float UE = towBack->get_UE(0).at(ieta);
0426 float tower_UE = UE * (1 + 2 * background_v2 * std::cos(2 * (tower_phi - background_Psi2)));
0427 tower_eT = (tower->get_energy() - tower_UE) / std::cosh(tower_eta);
0428 }
0429 }
0430
0431
0432 eFrac_emcal += tower_eT;
0433 jet_total_eT += tower_eT;
0434 }
0435 }
0436 struct jetty myjet;
0437 myjet.pt = jet->get_pt();
0438 myjet.eta = jet->get_eta();
0439 myjet.phi = jet->get_phi();
0440 myjet.et = jet_total_eT;
0441 myjet.emcal = eFrac_emcal;
0442 myjet.ihcal = eFrac_ihcal;
0443 myjet.ohcal = eFrac_ohcal;
0444 jets.push_back(myjet);
0445 }
0446 CLHEP::Hep3Vector vertex(0, 0, vtx_z);
0447 CLHEP::Hep3Vector max_cluster(0,0,0);
0448 std::vector<CLHEP::Hep3Vector> rawclusters{};
0449 if (clusters)
0450 {
0451
0452
0453 auto clusterrange = clusters->getClusters();
0454 for (auto iclus = clusterrange.first; iclus != clusterrange.second; ++iclus)
0455 {
0456
0457 RawCluster *cluster = iclus->second;
0458 CLHEP::Hep3Vector E_vec_cluster = RawClusterUtility::GetECoreVec(*cluster, vertex);
0459
0460 float energy = E_vec_cluster.mag();
0461
0462 if (energy < 1) continue;
0463 rawclusters.push_back(E_vec_cluster);
0464 if (energy > max_cluster.mag()) max_cluster = E_vec_cluster;
0465 }
0466 }
0467 std::string histoname = "_";
0468 bool isMB = false;
0469 if (m_useEmulator) isMB = true;
0470 else
0471 {
0472 isMB = triggeranalyzer->didTriggerFire("MBD N&S >= 1");
0473 }
0474 if (rawclusters.size() > 0)
0475 {
0476
0477 if (isMB)
0478 {
0479 auto h_cluster_max_et_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_max_et_mb", histoname.c_str())));
0480 h_cluster_max_et_mb->Fill(max_cluster.mag());
0481 auto h_cluster_max_et_scaled_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_max_et_scaled_mb", histoname.c_str())));
0482 h_cluster_max_et_scaled_mb->Fill(max_cluster.mag(), m_mbd_prescale);
0483 auto h_cluster_max_dist_mb = dynamic_cast<TH2*>(hm->getHisto(Form("h_2d%scluster_max_dist_mb", histoname.c_str())));
0484 h_cluster_max_dist_mb->Fill(max_cluster.pseudoRapidity(), max_cluster.phi());
0485 auto h_cluster_max_pt_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_max_pt_mb", histoname.c_str())));
0486 h_cluster_max_pt_mb->Fill(max_cluster.perp());
0487 auto h_cluster_max_pt_scaled_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_max_pt_scaled_mb", histoname.c_str())));
0488 h_cluster_max_pt_scaled_mb->Fill(max_cluster.perp(), m_mbd_prescale);
0489
0490 auto h_cluster_et_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_et_mb", histoname.c_str())));
0491 auto h_cluster_pt_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_pt_mb", histoname.c_str())));
0492 auto h_cluster_et_scaled_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_et_scaled_mb", histoname.c_str())));
0493 auto h_cluster_pt_scaled_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_pt_scaled_mb", histoname.c_str())));
0494 auto h_cluster_dist_mb = dynamic_cast<TH2*>(hm->getHisto(Form("h_2d%scluster_dist_mb", histoname.c_str())));
0495
0496 for (auto cluster : rawclusters)
0497 {
0498 h_cluster_et_mb->Fill(cluster.mag());
0499 h_cluster_et_scaled_mb->Fill(cluster.mag(), m_mbd_prescale);
0500 h_cluster_dist_mb->Fill(cluster.pseudoRapidity(), cluster.phi());
0501 h_cluster_pt_mb->Fill(cluster.perp());
0502 h_cluster_pt_scaled_mb->Fill(cluster.perp(), m_mbd_prescale);
0503 }
0504 for (int j = 0; j < 4; j++)
0505 {
0506 bool firerare = false;
0507
0508
0509 if (m_useEmulator)
0510 firerare = triggeranalyzer->didTriggerFire(m_photon_emu_triggernames[j]);
0511 else
0512 firerare = triggeranalyzer->checkRawTrigger(m_photon_triggernames[j]);
0513 if (firerare)
0514 {
0515 auto h_cluster_max_et_raw_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_max_et_raw_trig_%d", histoname.c_str(), j)));
0516 h_cluster_max_et_raw_trig->Fill(max_cluster.mag());
0517 auto h_cluster_max_pt_raw_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_max_pt_raw_trig_%d", histoname.c_str(), j)));
0518 h_cluster_max_pt_raw_trig->Fill(max_cluster.perp());
0519
0520 auto h_cluster_et_raw_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_et_raw_trig_%d", histoname.c_str(), j)));
0521 auto h_cluster_pt_raw_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_pt_raw_trig_%d", histoname.c_str(), j)));
0522 for (auto &cluster : rawclusters)
0523 {
0524 h_cluster_et_raw_trig->Fill(cluster.mag());
0525 h_cluster_pt_raw_trig->Fill(cluster.perp());
0526 }
0527
0528 }
0529
0530 }
0531 }
0532
0533 for (int j = 0; j < 4; j++)
0534 {
0535 if (m_useEmulator) continue;
0536 if (triggeranalyzer->didTriggerFire(m_photon_triggernames[j]))
0537 {
0538 if ( m_photon_prescales[j] == -1 ) m_photon_prescales[j] = triggeranalyzer->getTriggerPrescale(m_photon_triggernames[j]) + 1;
0539
0540 auto h_cluster_max_et_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_max_et_trig_%d", histoname.c_str(), j)));
0541 h_cluster_max_et_trig->Fill(max_cluster.mag());
0542 auto h_cluster_max_et_scaled_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_max_et_scaled_trig_%d", histoname.c_str(), j)));
0543 h_cluster_max_et_scaled_trig->Fill(max_cluster.mag(), m_photon_prescales[j]);
0544 auto h_cluster_max_dist_trig = dynamic_cast<TH2*>(hm->getHisto(Form("h_2d%scluster_max_dist_trig_%d", histoname.c_str(), j)));
0545 h_cluster_max_dist_trig->Fill(max_cluster.pseudoRapidity(), max_cluster.phi());
0546
0547 auto h_cluster_max_pt_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_max_pt_trig_%d", histoname.c_str(), j)));
0548 h_cluster_max_pt_trig->Fill(max_cluster.perp());
0549 auto h_cluster_max_pt_scaled_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_max_pt_scaled_trig_%d", histoname.c_str(), j)));
0550 h_cluster_max_pt_scaled_trig->Fill(max_cluster.perp(), m_photon_prescales[j]);
0551
0552
0553 auto h_cluster_et_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_et_trig_%d", histoname.c_str(), j)));
0554 auto h_cluster_pt_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_pt_trig_%d", histoname.c_str(), j)));
0555
0556 auto h_cluster_et_scaled_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_et_scaled_trig_%d", histoname.c_str(), j)));
0557 auto h_cluster_pt_scaled_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%scluster_pt_scaled_trig_%d", histoname.c_str(), j)));
0558
0559 auto h_cluster_dist_trig = dynamic_cast<TH2*>(hm->getHisto(Form("h_2d%scluster_dist_trig_%d", histoname.c_str(), j)));
0560
0561 for (auto &cluster : rawclusters)
0562 {
0563 h_cluster_et_trig->Fill(cluster.mag());
0564 h_cluster_et_scaled_trig->Fill(cluster.mag(), m_photon_prescales[j]);
0565 h_cluster_dist_trig->Fill(cluster.pseudoRapidity(), cluster.phi());
0566
0567 h_cluster_pt_trig->Fill(cluster.perp());
0568 h_cluster_pt_scaled_trig->Fill(cluster.perp(), m_photon_prescales[j]);
0569 }
0570 }
0571 }
0572 }
0573
0574 if (jets.size() == 0) return Fun4AllReturnCodes::EVENT_OK;
0575
0576 auto pt_ordered = jets;
0577 std::sort(jets.begin(), jets.end(), [](auto a, auto b) { return a.et > b.et; });
0578 std::sort(pt_ordered.begin(), pt_ordered.end(), [](auto a, auto b) { return a.pt > b.pt; });
0579
0580
0581 float max_et = jets.at(0).et;
0582 float max_pt = pt_ordered.at(0).pt;
0583 float max_eta = jets.at(0).eta;
0584 float max_phi = jets.at(0).phi;
0585
0586
0587 if (isMB)
0588 {
0589
0590 if (max_et > 4)
0591 {
0592 auto h_jet_max_et_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_max_et_mb", histoname.c_str())));
0593 h_jet_max_et_mb->Fill(max_et);
0594 auto h_jet_max_et_scaled_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_max_et_scaled_mb", histoname.c_str())));
0595 h_jet_max_et_scaled_mb->Fill(max_et, m_mbd_prescale);
0596
0597 auto h_jet_max_dist_mb = dynamic_cast<TH2*>(hm->getHisto(Form("h_2d%sjet_max_dist_mb", histoname.c_str())));
0598 h_jet_max_dist_mb->Fill(max_eta, max_phi);
0599
0600 }
0601 if (max_pt > 4)
0602 {
0603 auto h_jet_max_pt_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_max_pt_mb", histoname.c_str())));
0604 h_jet_max_pt_mb->Fill(max_pt);
0605 auto h_jet_max_pt_scaled_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_max_pt_scaled_mb", histoname.c_str())));
0606 h_jet_max_pt_scaled_mb->Fill(max_pt, m_mbd_prescale);
0607
0608 }
0609
0610 auto h_jet_et_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_et_mb", histoname.c_str())));
0611 auto h_jet_pt_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_pt_mb", histoname.c_str())));
0612 auto h_jet_et_scaled_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_et_scaled_mb", histoname.c_str())));
0613 auto h_jet_pt_scaled_mb = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_pt_scaled_mb", histoname.c_str())));
0614 auto h_jet_dist_mb = dynamic_cast<TH2*>(hm->getHisto(Form("h_2d%sjet_dist_mb", histoname.c_str())));
0615
0616 for (auto &jet : jets)
0617 {
0618 if (max_et > 4)
0619 {
0620 h_jet_et_mb->Fill(jet.et);
0621 h_jet_et_scaled_mb->Fill(jet.et, m_mbd_prescale);
0622 h_jet_dist_mb->Fill(jet.eta, jet.phi);
0623 }
0624 if (max_et > 4)
0625 {
0626 h_jet_pt_mb->Fill(jet.pt);
0627 h_jet_pt_scaled_mb->Fill(jet.pt, m_mbd_prescale);
0628 }
0629 }
0630 for (int j = 0; j < 4; j++)
0631 {
0632 bool firerare = false;
0633
0634
0635 if (m_useEmulator)
0636 firerare = triggeranalyzer->didTriggerFire(m_jet_emu_triggernames[j]);
0637 else
0638 firerare = triggeranalyzer->checkRawTrigger(m_jet_triggernames[j]);
0639 if (firerare)
0640 {
0641 if (max_et > 4)
0642 {
0643 auto h_jet_max_et_raw_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_max_et_raw_trig_%d", histoname.c_str(), j)));
0644 h_jet_max_et_raw_trig->Fill(max_et);
0645 }
0646 if (max_pt > 4)
0647 {
0648
0649 auto h_jet_max_pt_raw_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_max_pt_raw_trig_%d", histoname.c_str(), j)));
0650 h_jet_max_pt_raw_trig->Fill(max_pt);
0651 }
0652
0653 auto h_jet_et_raw_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_et_raw_trig_%d", histoname.c_str(), j)));
0654 auto h_jet_pt_raw_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_pt_raw_trig_%d", histoname.c_str(), j)));
0655 for (auto &jet : jets)
0656 {
0657 if (max_et > 4)
0658 h_jet_et_raw_trig->Fill(jet.et);
0659 if (max_pt > 4)
0660 h_jet_pt_raw_trig->Fill(jet.pt);
0661 }
0662
0663 }
0664
0665 }
0666 }
0667
0668 for (int j = 0; j < 4; j++)
0669 {
0670 if (m_useEmulator) continue;
0671 if (triggeranalyzer->didTriggerFire(m_jet_triggernames[j]))
0672 {
0673 if ( m_jet_prescales[j] == -1 ) m_jet_prescales[j] = triggeranalyzer->getTriggerPrescale(m_jet_triggernames[j]) + 1;
0674
0675 if (max_et > 4)
0676 {
0677 auto h_jet_max_et_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_max_et_trig_%d", histoname.c_str(), j)));
0678 h_jet_max_et_trig->Fill(max_et);
0679 auto h_jet_max_et_scaled_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_max_et_scaled_trig_%d", histoname.c_str(), j)));
0680 h_jet_max_et_scaled_trig->Fill(max_et, m_jet_prescales[j]);
0681 auto h_jet_max_dist_trig = dynamic_cast<TH2*>(hm->getHisto(Form("h_2d%sjet_max_dist_trig_%d", histoname.c_str(), j)));
0682 h_jet_max_dist_trig->Fill(max_eta, max_phi);
0683 }
0684 if (max_pt > 4)
0685 {
0686 auto h_jet_max_pt_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_max_pt_trig_%d", histoname.c_str(), j)));
0687 h_jet_max_pt_trig->Fill(max_pt);
0688 auto h_jet_max_pt_scaled_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_max_pt_scaled_trig_%d", histoname.c_str(), j)));
0689 h_jet_max_pt_scaled_trig->Fill(max_pt, m_jet_prescales[j]);
0690 }
0691
0692 auto h_jet_et_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_et_trig_%d", histoname.c_str(), j)));
0693 auto h_jet_pt_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_pt_trig_%d", histoname.c_str(), j)));
0694
0695 auto h_jet_et_scaled_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_et_scaled_trig_%d", histoname.c_str(), j)));
0696 auto h_jet_pt_scaled_trig = dynamic_cast<TH1*>(hm->getHisto(Form("h%sjet_pt_scaled_trig_%d", histoname.c_str(), j)));
0697
0698 auto h_jet_dist_trig = dynamic_cast<TH2*>(hm->getHisto(Form("h_2d%sjet_dist_trig_%d", histoname.c_str(), j)));
0699
0700 for (auto &jet : jets)
0701 {
0702 if (max_et > 4)
0703 {
0704 h_jet_et_trig->Fill(jet.et);
0705 h_jet_et_scaled_trig->Fill(jet.et, m_jet_prescales[j]);
0706 h_jet_dist_trig->Fill(jet.eta, jet.phi);
0707 }
0708 if (max_pt > 4)
0709 {
0710 h_jet_pt_trig->Fill(jet.pt);
0711 h_jet_pt_scaled_trig->Fill(jet.pt, m_jet_prescales[j]);
0712 }
0713 }
0714 }
0715 }
0716
0717
0718 return Fun4AllReturnCodes::EVENT_OK;
0719 }
0720
0721
0722
0723 void TriggerPerformance::GetNodes(PHCompositeNode* topNode)
0724 {
0725
0726
0727 }
0728
0729 int TriggerPerformance::ResetEvent(PHCompositeNode *topNode)
0730 {
0731
0732 return Fun4AllReturnCodes::EVENT_OK;
0733 }
0734
0735
0736 int TriggerPerformance::EndRun(const int runnumber)
0737 {
0738 return Fun4AllReturnCodes::EVENT_OK;
0739 }
0740
0741
0742 int TriggerPerformance::End(PHCompositeNode *topNode)
0743 {
0744 if (Verbosity() > 0)
0745 {
0746 std::cout << "TriggerPerformance::End(PHCompositeNode *topNode) This is the End..." << std::endl;
0747 }
0748 hm->dumpHistos(_foutname.c_str());
0749 return Fun4AllReturnCodes::EVENT_OK;
0750 }
0751
0752
0753 int TriggerPerformance::Reset(PHCompositeNode *topNode)
0754 {
0755 if (Verbosity() > 0)
0756 {
0757 std::cout << "TriggerPerformance::Reset(PHCompositeNode *topNode) being Reset" << std::endl;
0758 }
0759 return Fun4AllReturnCodes::EVENT_OK;
0760 }