Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:17:54

0001 #include "CaloEvaluator.h"
0002 
0003 #include "CaloEvalStack.h"
0004 #include "CaloRawClusterEval.h"
0005 #include "CaloRawTowerEval.h"
0006 #include "CaloTruthEval.h"
0007 
0008 #include <g4main/PHG4Particle.h>
0009 #include <g4main/PHG4Shower.h>
0010 #include <g4main/PHG4TruthInfoContainer.h>
0011 #include <g4main/PHG4VtxPoint.h>
0012 
0013 #include <globalvertex/GlobalVertex.h>
0014 #include <globalvertex/GlobalVertexMap.h>
0015 
0016 #include <calobase/RawCluster.h>
0017 #include <calobase/RawClusterContainer.h>
0018 #include <calobase/RawClusterUtility.h>
0019 #include <calobase/RawTower.h>
0020 #include <calobase/RawTowerContainer.h>
0021 #include <calobase/RawTowerGeom.h>
0022 #include <calobase/RawTowerGeomContainer.h>
0023 #include <calobase/TowerInfo.h>
0024 #include <calobase/TowerInfoContainer.h>
0025 
0026 #include <fun4all/Fun4AllReturnCodes.h>
0027 #include <fun4all/SubsysReco.h>
0028 
0029 #include <phool/getClass.h>
0030 #include <phool/phool.h>
0031 
0032 #include <TFile.h>
0033 #include <TNtuple.h>
0034 #include <TTree.h>
0035 
0036 #include <CLHEP/Vector/ThreeVector.h>
0037 
0038 #include <cmath>
0039 #include <cstdlib>
0040 #include <iostream>
0041 #include <set>
0042 #include <utility>
0043 
0044 CaloEvaluator::CaloEvaluator(const std::string& name, const std::string& caloname, const std::string& filename)
0045   : SubsysReco(name)
0046   , _caloname(caloname)
0047   , _filename(filename)
0048 {
0049 }
0050 
0051 int CaloEvaluator::Init(PHCompositeNode* /*topNode*/)
0052 {
0053   delete _tfile;  // make cppcheck happy
0054   _tfile = new TFile(_filename.c_str(), "RECREATE");
0055 
0056   if (_do_gpoint_eval)
0057   {
0058     _ntp_gpoint = new TNtuple("ntp_gpoint", "primary vertex => best (first) vertex",
0059                               "event:gvx:gvy:gvz:"
0060                               "vx:vy:vz");
0061   }
0062 
0063   if (_do_gshower_eval)
0064   {
0065     _ntp_gshower = new TNtuple("ntp_gshower", "truth shower => best cluster",
0066                                "event:gparticleID:gflavor:gnhits:"
0067                                "geta:gphi:ge:gpt:gvx:gvy:gvz:gembed:gedep:"
0068                                "clusterID:ntowers:eta:x:y:z:phi:e:efromtruth");
0069   }
0070 
0071   // Barak: Added TTree to will allow the TowerID to be set correctly as integer
0072   if (_do_tower_eval)
0073   {
0074     _ntp_tower = new TNtuple("ntp_tower", "tower => max truth primary",
0075                              "event:towerID:ieta:iphi:eta:phi:e:x:y:z:"
0076                              "gparticleID:gflavor:gnhits:"
0077                              "geta:gphi:ge:gpt:gvx:gvy:gvz:"
0078                              "gembed:gedep:"
0079                              "efromtruth");
0080 
0081     // Make Tree
0082     _tower_debug = new TTree("tower_debug", "tower => max truth primary");
0083 
0084     _tower_debug->Branch("event", &_ievent, "event/I");
0085     _tower_debug->Branch("towerID", &_towerID_debug, "towerID/I");
0086     _tower_debug->Branch("ieta", &_ieta_debug, "ieta/I");
0087     _tower_debug->Branch("iphi", &_iphi_debug, "iphi/I");
0088     _tower_debug->Branch("eta", &_eta_debug, "eta/F");
0089     _tower_debug->Branch("phi", &_phi_debug, "phi/F");
0090     _tower_debug->Branch("e", &_e_debug, "e/F");
0091     _tower_debug->Branch("x", &_x_debug, "x/F");
0092     _tower_debug->Branch("y", &_y_debug, "y/F");
0093     _tower_debug->Branch("z", &_z_debug, "z/F");
0094   }
0095 
0096   if (_do_cluster_eval)
0097   {
0098     _ntp_cluster = new TNtuple("ntp_cluster", "cluster => max truth primary",
0099                                "event:clusterID:ntowers:eta:x:y:z:phi:e:"
0100                                "gparticleID:gflavor:gnhits:"
0101                                "geta:gphi:ge:gpt:gvx:gvy:gvz:"
0102                                "gembed:gedep:"
0103                                "efromtruth");
0104   }
0105 
0106   return Fun4AllReturnCodes::EVENT_OK;
0107 }
0108 
0109 int CaloEvaluator::process_event(PHCompositeNode* topNode)
0110 {
0111   if (!_caloevalstack)
0112   {
0113     _caloevalstack = new CaloEvalStack(topNode, _caloname);
0114     _caloevalstack->set_strict(_strict);
0115     _caloevalstack->set_verbosity(Verbosity() + 1);
0116   }
0117   else
0118   {
0119     _caloevalstack->next_event(topNode);
0120   }
0121 
0122   //-----------------------------------
0123   // print what is coming into the code
0124   //-----------------------------------
0125 
0126   printInputInfo(topNode);
0127 
0128   //---------------------------
0129   // fill the Evaluator NTuples
0130   //---------------------------
0131 
0132   fillOutputNtuples(topNode);
0133 
0134   //--------------------------------------------------
0135   // Print out the ancestry information for this event
0136   //--------------------------------------------------
0137 
0138   printOutputInfo(topNode);
0139 
0140   ++_ievent;
0141   ++m_EvtCounter;
0142 
0143   return Fun4AllReturnCodes::EVENT_OK;
0144 }
0145 
0146 int CaloEvaluator::End(PHCompositeNode* /*topNode*/)
0147 {
0148   _tfile->cd();
0149 
0150   if (_do_gpoint_eval)
0151   {
0152     _ntp_gpoint->Write();
0153   }
0154   if (_do_gshower_eval)
0155   {
0156     _ntp_gshower->Write();
0157   }
0158   if (_do_tower_eval)
0159   {
0160     _ntp_tower->Write();
0161     _tower_debug->Write();
0162   }
0163   if (_do_cluster_eval)
0164   {
0165     _ntp_cluster->Write();
0166   }
0167 
0168   _tfile->Close();
0169 
0170   delete _tfile;
0171 
0172   if (Verbosity() > 0)
0173   {
0174     std::cout << "========================= " << Name() << "::End() ============================" << std::endl;
0175     std::cout << " " << m_EvtCounter << " events of output written to: " << _filename << std::endl;
0176     std::cout << "===========================================================================" << std::endl;
0177   }
0178 
0179   if (_caloevalstack)
0180   {
0181     delete _caloevalstack;
0182   }
0183 
0184   return Fun4AllReturnCodes::EVENT_OK;
0185 }
0186 
0187 void CaloEvaluator::printInputInfo(PHCompositeNode* topNode)
0188 {
0189   if (Verbosity() > 2)
0190   {
0191     std::cout << "CaloEvaluator::printInputInfo() entered" << std::endl;
0192   }
0193 
0194   // print out the truth container
0195 
0196   if (Verbosity() > 1)
0197   {
0198     std::cout << std::endl;
0199     std::cout << PHWHERE << "   NEW INPUT FOR EVENT " << _ievent << std::endl;
0200     std::cout << std::endl;
0201 
0202     // need things off of the DST...
0203     PHG4TruthInfoContainer* truthinfo = findNode::getClass<PHG4TruthInfoContainer>(topNode, "G4TruthInfo");
0204     if (!truthinfo)
0205     {
0206       std::cout << PHWHERE << " ERROR: Can't find G4TruthInfo" << std::endl;
0207       exit(-1);
0208     }
0209 
0210     std::cout << Name() << ": PHG4TruthInfoContainer contents: " << std::endl;
0211 
0212     PHG4TruthInfoContainer::Range truthrange = truthinfo->GetParticleRange();
0213     for (PHG4TruthInfoContainer::Iterator truthiter = truthrange.first;
0214          truthiter != truthrange.second;
0215          ++truthiter)
0216     {
0217       PHG4Particle* particle = truthiter->second;
0218 
0219       std::cout << truthiter->first << " => pid: " << particle->get_pid()
0220                 << " pt: " << sqrt(pow(particle->get_px(), 2) + pow(particle->get_py(), 2)) << std::endl;
0221     }
0222   }
0223 
0224   return;
0225 }
0226 
0227 void CaloEvaluator::printOutputInfo(PHCompositeNode* topNode)
0228 {
0229   if (Verbosity() > 2)
0230   {
0231     std::cout << "CaloEvaluator::printOutputInfo() entered" << std::endl;
0232   }
0233 
0234   CaloRawClusterEval* clustereval = _caloevalstack->get_rawcluster_eval();
0235   clustereval->set_usetowerinfo(_use_towerinfo);
0236   CaloTruthEval* trutheval = _caloevalstack->get_truth_eval();
0237 
0238   //==========================================
0239   // print out some useful stuff for debugging
0240   //==========================================
0241 
0242   if (Verbosity() > 1)
0243   {
0244     // event information
0245     std::cout << std::endl;
0246     std::cout << PHWHERE << "   NEW OUTPUT FOR EVENT " << _ievent << std::endl;
0247     std::cout << std::endl;
0248 
0249     // need things off of the DST...
0250     PHG4TruthInfoContainer* truthinfo = findNode::getClass<PHG4TruthInfoContainer>(topNode, "G4TruthInfo");
0251     if (!truthinfo)
0252     {
0253       std::cout << PHWHERE << " ERROR: Can't find G4TruthInfo" << std::endl;
0254       exit(-1);
0255     }
0256 
0257     // need things off of the DST...
0258     GlobalVertexMap* vertexmap = findNode::getClass<GlobalVertexMap>(topNode, "GlobalVertexMap");
0259 
0260     PHG4VtxPoint* gvertex = truthinfo->GetPrimaryVtx(truthinfo->GetPrimaryVertexIndex());
0261     float gvx = gvertex->get_x();
0262     float gvy = gvertex->get_y();
0263     float gvz = gvertex->get_z();
0264 
0265     float vx = NAN;
0266     float vy = NAN;
0267     float vz = NAN;
0268     if (vertexmap)
0269     {
0270       if (!vertexmap->empty())
0271       {
0272         GlobalVertex* vertex = (vertexmap->begin()->second);
0273 
0274         vx = vertex->get_x();
0275         vy = vertex->get_y();
0276         vz = vertex->get_z();
0277       }
0278     }
0279 
0280     std::cout << "vtrue = (" << gvx << "," << gvy << "," << gvz << ") => vreco = (" << vx << "," << vy << "," << vz << ")" << std::endl;
0281 
0282     PHG4TruthInfoContainer::ConstRange range = truthinfo->GetPrimaryParticleRange();
0283     for (PHG4TruthInfoContainer::ConstIterator iter = range.first;
0284          iter != range.second;
0285          ++iter)
0286     {
0287       PHG4Particle* primary = iter->second;
0288 
0289       std::cout << std::endl;
0290 
0291       std::cout << "===Primary PHG4Particle=========================================" << std::endl;
0292       std::cout << " particle id = " << primary->get_track_id() << std::endl;
0293       std::cout << " flavor = " << primary->get_pid() << std::endl;
0294       std::cout << " (px,py,pz,e) = (";
0295 
0296       float gpx = primary->get_px();
0297       float gpy = primary->get_py();
0298       float gpz = primary->get_pz();
0299       float ge = primary->get_e();
0300 
0301       std::cout.width(5);
0302       std::cout << gpx;
0303       std::cout << ",";
0304       std::cout.width(5);
0305       std::cout << gpy;
0306       std::cout << ",";
0307       std::cout.width(5);
0308       std::cout << gpz;
0309       std::cout << ",";
0310       std::cout.width(5);
0311       std::cout << ge;
0312       std::cout << ")" << std::endl;
0313 
0314       float gpt = std::sqrt(gpx * gpx + gpy * gpy);
0315       float geta = NAN;
0316       if (gpt != 0.0)
0317       {
0318         geta = std::asinh(gpz / gpt);
0319       }
0320       float gphi = std::atan2(gpy, gpx);
0321 
0322       std::cout << "(eta,phi,e,pt) = (";
0323       std::cout.width(5);
0324       std::cout << geta;
0325       std::cout << ",";
0326       std::cout.width(5);
0327       std::cout << gphi;
0328       std::cout << ",";
0329       std::cout.width(5);
0330       std::cout << ge;
0331       std::cout << ",";
0332       std::cout.width(5);
0333       std::cout << gpt;
0334       std::cout << ")" << std::endl;
0335 
0336       PHG4VtxPoint* vtx = trutheval->get_vertex(primary);
0337       float local_gvx = vtx->get_x();
0338       float local_gvy = vtx->get_y();
0339       float local_gvz = vtx->get_z();
0340 
0341       std::cout << " vtrue = (";
0342       std::cout.width(5);
0343       std::cout << local_gvx;
0344       std::cout << ",";
0345       std::cout.width(5);
0346       std::cout << local_gvy;
0347       std::cout << ",";
0348       std::cout.width(5);
0349       std::cout << local_gvz;
0350       std::cout << ")" << std::endl;
0351 
0352       std::cout << " embed = " << trutheval->get_embed(primary) << std::endl;
0353       std::cout << " edep = " << trutheval->get_shower_energy_deposit(primary) << std::endl;
0354 
0355       std::set<RawCluster*> clusters = clustereval->all_clusters_from(primary);
0356       for (auto cluster : clusters)
0357       {
0358         float ntowers = cluster->getNTowers();
0359         float x = cluster->get_x();
0360         float y = cluster->get_y();
0361         float z = cluster->get_z();
0362         float phi = cluster->get_phi();
0363         float e = cluster->get_energy();
0364 
0365         float efromtruth = clustereval->get_energy_contribution(cluster, primary);
0366 
0367         std::cout << " => #" << cluster->get_id() << " (x,y,z,phi,e) = (";
0368         std::cout.width(5);
0369         std::cout << x;
0370         std::cout << ",";
0371         std::cout.width(5);
0372         std::cout << y;
0373         std::cout << ",";
0374         std::cout.width(5);
0375         std::cout << z;
0376         std::cout << ",";
0377         std::cout.width(5);
0378         std::cout << phi;
0379         std::cout << ",";
0380         std::cout.width(5);
0381         std::cout << e;
0382         std::cout << "), ntowers = " << ntowers << ", efromtruth = " << efromtruth << std::endl;
0383       }
0384     }
0385     std::cout << std::endl;
0386   }
0387 
0388   return;
0389 }
0390 
0391 void CaloEvaluator::fillOutputNtuples(PHCompositeNode* topNode)
0392 {
0393   if (Verbosity() > 2)
0394   {
0395     std::cout << "CaloEvaluator::fillOutputNtuples() entered" << std::endl;
0396   }
0397 
0398   CaloRawClusterEval* clustereval = _caloevalstack->get_rawcluster_eval();
0399   CaloRawTowerEval* towereval = _caloevalstack->get_rawtower_eval();
0400   CaloTruthEval* trutheval = _caloevalstack->get_truth_eval();
0401 
0402   //----------------------
0403   // fill the Event NTuple
0404   //----------------------
0405 
0406   if (_do_gpoint_eval)
0407   {
0408     // need things off of the DST...
0409     PHG4TruthInfoContainer* truthinfo = findNode::getClass<PHG4TruthInfoContainer>(topNode, "G4TruthInfo");
0410     if (!truthinfo)
0411     {
0412       std::cout << PHWHERE << " ERROR: Can't find G4TruthInfo" << std::endl;
0413       exit(-1);
0414     }
0415 
0416     // need things off of the DST...
0417     GlobalVertexMap* vertexmap = findNode::getClass<GlobalVertexMap>(topNode, "GlobalVertexMap");
0418 
0419     PHG4VtxPoint* gvertex = truthinfo->GetPrimaryVtx(truthinfo->GetPrimaryVertexIndex());
0420     float gvx = gvertex->get_x();
0421     float gvy = gvertex->get_y();
0422     float gvz = gvertex->get_z();
0423 
0424     float vx = NAN;
0425     float vy = NAN;
0426     float vz = NAN;
0427     if (vertexmap)
0428     {
0429       if (!vertexmap->empty())
0430       {
0431         GlobalVertex* vertex = (vertexmap->begin()->second);
0432 
0433         vx = vertex->get_x();
0434         vy = vertex->get_y();
0435         vz = vertex->get_z();
0436       }
0437     }
0438 
0439     float gpoint_data[7] = {(float) _ievent,
0440                             gvx,
0441                             gvy,
0442                             gvz,
0443                             vx,
0444                             vy,
0445                             vz};
0446 
0447     _ntp_gpoint->Fill(gpoint_data);
0448   }
0449 
0450   //------------------------
0451   // fill the Gshower NTuple
0452   //------------------------
0453 
0454   if (_ntp_gshower)
0455   {
0456     if (Verbosity() > 1)
0457     {
0458       std::cout << Name() << " CaloEvaluator::filling gshower ntuple..." << std::endl;
0459     }
0460 
0461     GlobalVertexMap* vertexmap = findNode::getClass<GlobalVertexMap>(topNode, "GlobalVertexMap");
0462 
0463     PHG4TruthInfoContainer* truthinfo = findNode::getClass<PHG4TruthInfoContainer>(topNode, "G4TruthInfo");
0464     if (!truthinfo)
0465     {
0466       std::cout << PHWHERE << " ERROR: Can't find G4TruthInfo" << std::endl;
0467       exit(-1);
0468     }
0469 
0470     PHG4TruthInfoContainer::ConstRange range = truthinfo->GetPrimaryParticleRange();
0471     for (PHG4TruthInfoContainer::ConstIterator iter = range.first;
0472          iter != range.second;
0473          ++iter)
0474     {
0475       PHG4Particle* primary = iter->second;
0476 
0477       if (primary->get_e() < _truth_e_threshold)
0478       {
0479         continue;
0480       }
0481 
0482       if (!_truth_trace_embed_flags.empty())
0483       {
0484         if (_truth_trace_embed_flags.find(trutheval->get_embed(primary)) ==
0485             _truth_trace_embed_flags.end())
0486         {
0487           continue;
0488         }
0489       }
0490 
0491       float gparticleID = primary->get_track_id();
0492       float gflavor = primary->get_pid();
0493 
0494       PHG4Shower* shower = trutheval->get_primary_shower(primary);
0495       float gnhits = NAN;
0496       if (shower)
0497       {
0498         gnhits = shower->get_nhits(trutheval->get_caloid());
0499       }
0500       else
0501       {
0502         gnhits = 0.0;
0503       }
0504       float gpx = primary->get_px();
0505       float gpy = primary->get_py();
0506       float gpz = primary->get_pz();
0507       float ge = primary->get_e();
0508 
0509       float gpt = std::sqrt(gpx * gpx + gpy * gpy);
0510       float geta = NAN;
0511       if (gpt != 0.0)
0512       {
0513         geta = std::asinh(gpz / gpt);
0514       }
0515       float gphi = std::atan2(gpy, gpx);
0516 
0517       PHG4VtxPoint* vtx = trutheval->get_vertex(primary);
0518       float gvx = vtx->get_x();
0519       float gvy = vtx->get_y();
0520       float gvz = vtx->get_z();
0521 
0522       float gembed = trutheval->get_embed(primary);
0523       float gedep = trutheval->get_shower_energy_deposit(primary);
0524 
0525       RawCluster* cluster = clustereval->best_cluster_from(primary);
0526 
0527       float clusterID = NAN;
0528       float ntowers = NAN;
0529       float eta = NAN;
0530       float x = NAN;
0531       float y = NAN;
0532       float z = NAN;
0533       float phi = NAN;
0534       float e = NAN;
0535 
0536       float efromtruth = NAN;
0537 
0538       if (cluster)
0539       {
0540         clusterID = cluster->get_id();
0541         ntowers = cluster->getNTowers();
0542         x = cluster->get_x();
0543         y = cluster->get_y();
0544         z = cluster->get_z();
0545         phi = cluster->get_phi();
0546         e = cluster->get_energy();
0547 
0548         efromtruth = clustereval->get_energy_contribution(cluster, primary);
0549 
0550         // require vertex for cluster eta calculation
0551         if (vertexmap)
0552         {
0553           if (!vertexmap->empty())
0554           {
0555             GlobalVertex* vertex = (vertexmap->begin()->second);
0556 
0557             eta =
0558                 RawClusterUtility::GetPseudorapidity(
0559                     *cluster,
0560                     CLHEP::Hep3Vector(vertex->get_x(), vertex->get_y(), vertex->get_z()));
0561           }
0562         }
0563       }
0564 
0565       float shower_data[] = {(float) _ievent,
0566                              gparticleID,
0567                              gflavor,
0568                              gnhits,
0569                              geta,
0570                              gphi,
0571                              ge,
0572                              gpt,
0573                              gvx,
0574                              gvy,
0575                              gvz,
0576                              gembed,
0577                              gedep,
0578                              clusterID,
0579                              ntowers,
0580                              eta,
0581                              x,
0582                              y,
0583                              z,
0584                              phi,
0585                              e,
0586                              efromtruth};
0587 
0588       _ntp_gshower->Fill(shower_data);
0589     }
0590   }
0591 
0592   //----------------------
0593   // fill the Tower NTuple
0594   //----------------------
0595 
0596   if (_do_tower_eval)
0597   {
0598     if (!_use_towerinfo)
0599     {
0600       if (Verbosity() > 1)
0601       {
0602         std::cout << "CaloEvaluator::filling tower ntuple..." << std::endl;
0603       }
0604 
0605       std::string towernode = "TOWER_CALIB_" + _caloname;
0606       RawTowerContainer* towers = findNode::getClass<RawTowerContainer>(topNode, towernode.c_str());
0607       if (!towers)
0608       {
0609         std::cout << PHWHERE << " ERROR: Can't find " << towernode << std::endl;
0610         exit(-1);
0611       }
0612 
0613       std::string towergeomnode = "TOWERGEOM_" + _caloname;
0614       RawTowerGeomContainer* towergeom = findNode::getClass<RawTowerGeomContainer>(topNode, towergeomnode.c_str());
0615       if (!towergeom)
0616       {
0617         std::cout << PHWHERE << " ERROR: Can't find " << towergeomnode << std::endl;
0618         exit(-1);
0619       }
0620 
0621       RawTowerContainer::ConstRange begin_end = towers->getTowers();
0622       RawTowerContainer::ConstIterator rtiter;
0623       for (rtiter = begin_end.first; rtiter != begin_end.second; ++rtiter)
0624       {
0625         RawTower* tower = rtiter->second;
0626 
0627         if (tower->get_energy() < _reco_e_threshold)
0628         {
0629           continue;
0630         }
0631 
0632         RawTowerGeom* tower_geom = towergeom->get_tower_geometry(tower->get_id());
0633         if (!tower_geom)
0634         {
0635           std::cout << PHWHERE << " ERROR: Can't find tower geometry for this tower hit: ";
0636           tower->identify();
0637           exit(-1);
0638         }
0639 
0640         // std::cout<<"Tower ID = "<<tower->get_id()<<" for bin(j,k)= "<<tower->get_bineta()<<","<<tower->get_binphi()<<std::endl; //Added by Barak
0641         const float towerid = tower->get_id();
0642         const float ieta = tower->get_bineta();
0643         const float iphi = tower->get_binphi();
0644         const float eta = tower_geom->get_eta();
0645         const float phi = tower_geom->get_phi();
0646         const float e = tower->get_energy();
0647         const float x = tower_geom->get_center_x();
0648         const float y = tower_geom->get_center_y();
0649         const float z = tower_geom->get_center_z();
0650 
0651         // Added by Barak
0652         _towerID_debug = tower->get_id();
0653         _ieta_debug = tower->get_bineta();
0654         _iphi_debug = tower->get_binphi();
0655         _eta_debug = tower_geom->get_eta();
0656         _phi_debug = tower_geom->get_phi();
0657         _e_debug = tower->get_energy();
0658         _x_debug = tower_geom->get_center_x();
0659         _y_debug = tower_geom->get_center_y();
0660         _z_debug = tower_geom->get_center_z();
0661 
0662         PHG4Particle* primary = towereval->max_truth_primary_particle_by_energy(tower);
0663 
0664         float gparticleID = NAN;
0665         float gflavor = NAN;
0666         float gnhits = NAN;
0667         float gpx = NAN;
0668         float gpy = NAN;
0669         float gpz = NAN;
0670         float ge = NAN;
0671 
0672         float gpt = NAN;
0673         float geta = NAN;
0674         float gphi = NAN;
0675 
0676         float gvx = NAN;
0677         float gvy = NAN;
0678         float gvz = NAN;
0679 
0680         float gembed = NAN;
0681         float gedep = NAN;
0682 
0683         float efromtruth = NAN;
0684 
0685         if (primary)
0686         {
0687           gparticleID = primary->get_track_id();
0688           gflavor = primary->get_pid();
0689 
0690           PHG4Shower* shower = trutheval->get_primary_shower(primary);
0691           if (shower)
0692           {
0693             gnhits = shower->get_nhits(trutheval->get_caloid());
0694           }
0695           else
0696           {
0697             gnhits = 0.0;
0698           }
0699           gpx = primary->get_px();
0700           gpy = primary->get_py();
0701           gpz = primary->get_pz();
0702           ge = primary->get_e();
0703 
0704           gpt = std::sqrt(gpx * gpx + gpy * gpy);
0705           if (gpt != 0.0)
0706           {
0707             geta = std::asinh(gpz / gpt);
0708           }
0709           gphi = std::atan2(gpy, gpx);
0710 
0711           PHG4VtxPoint* vtx = trutheval->get_vertex(primary);
0712 
0713           if (vtx)
0714           {
0715             gvx = vtx->get_x();
0716             gvy = vtx->get_y();
0717             gvz = vtx->get_z();
0718           }
0719 
0720           gembed = trutheval->get_embed(primary);
0721           gedep = trutheval->get_shower_energy_deposit(primary);
0722 
0723           efromtruth = towereval->get_energy_contribution(tower, primary);
0724         }
0725 
0726         float tower_data[] = {(float) _ievent,
0727                               towerid,
0728                               ieta,
0729                               iphi,
0730                               eta,
0731                               phi,
0732                               e,
0733                               x,
0734                               y,
0735                               z,
0736                               gparticleID,
0737                               gflavor,
0738                               gnhits,
0739                               geta,
0740                               gphi,
0741                               ge,
0742                               gpt,
0743                               gvx,
0744                               gvy,
0745                               gvz,
0746                               gembed,
0747                               gedep,
0748                               efromtruth};
0749 
0750         _ntp_tower->Fill(tower_data);
0751         _tower_debug->Fill();  // Added by Barak (see above for explanation)
0752       }
0753     }
0754     else
0755     {
0756       if (Verbosity() > 1)
0757       {
0758         std::cout << "CaloEvaluator::filling tower ntuple..." << std::endl;
0759       }
0760 
0761       std::string towernode = "TOWERINFO_CALIB_" + _caloname;
0762       TowerInfoContainer* towers = findNode::getClass<TowerInfoContainer>(topNode, towernode.c_str());
0763       if (!towers)
0764       {
0765         std::cout << PHWHERE << " ERROR: Can't find " << towernode << std::endl;
0766         exit(-1);
0767       }
0768 
0769       std::string towergeomnode = "TOWERGEOM_" + _caloname;
0770       RawTowerGeomContainer* towergeom = findNode::getClass<RawTowerGeomContainer>(topNode, towergeomnode.c_str());
0771       if (!towergeom)
0772       {
0773         std::cout << PHWHERE << " ERROR: Can't find " << towergeomnode << std::endl;
0774         exit(-1);
0775       }
0776 
0777       unsigned int ntowers = towers->size();
0778 
0779       for (unsigned int channel = 0; channel < ntowers; channel++)
0780       {
0781         TowerInfo* tower = towers->get_tower_at_channel(channel);
0782 
0783         if (tower->get_energy() < _reco_e_threshold)
0784         {
0785           continue;
0786         }
0787 
0788         unsigned int towerkey = towers->encode_key(channel);
0789         unsigned int ieta_bin = towers->getTowerEtaBin(towerkey);
0790         unsigned int iphi_bin = towers->getTowerPhiBin(towerkey);
0791         const RawTowerDefs::keytype key = RawTowerDefs::encode_towerid(
0792             RawTowerDefs::convert_name_to_caloid("CEMC"), ieta_bin, iphi_bin);
0793         RawTowerGeom* tower_geom = towergeom->get_tower_geometry(key);
0794         if (!tower_geom)
0795         {
0796           std::cout << PHWHERE << " ERROR: Can't find tower geometry for this tower hit: ";
0797           tower->identify();
0798           exit(-1);
0799         }
0800 
0801         // std::cout<<"Tower ID = "<<tower->get_id()<<" for bin(j,k)= "<<tower->get_bineta()<<","<<tower->get_binphi()<<std::endl; //Added by Barak
0802         const float towerid = key;
0803         const float ieta = ieta_bin;
0804         const float iphi = iphi_bin;
0805         const float eta = tower_geom->get_eta();
0806         const float phi = tower_geom->get_phi();
0807         const float e = tower->get_energy();
0808         const float x = tower_geom->get_center_x();
0809         const float y = tower_geom->get_center_y();
0810         const float z = tower_geom->get_center_z();
0811 
0812         // Added by Barak
0813         _towerID_debug = towerid;
0814         _ieta_debug = ieta;
0815         _iphi_debug = iphi;
0816         _eta_debug = eta;
0817         _phi_debug = phi;
0818         _e_debug = e;
0819         _x_debug = x;
0820         _y_debug = y;
0821         _z_debug = z;
0822 
0823         PHG4Particle* primary = towereval->max_truth_primary_particle_by_energy(tower);
0824 
0825         float gparticleID = NAN;
0826         float gflavor = NAN;
0827         float gnhits = NAN;
0828         float gpx = NAN;
0829         float gpy = NAN;
0830         float gpz = NAN;
0831         float ge = NAN;
0832 
0833         float gpt = NAN;
0834         float geta = NAN;
0835         float gphi = NAN;
0836 
0837         float gvx = NAN;
0838         float gvy = NAN;
0839         float gvz = NAN;
0840 
0841         float gembed = NAN;
0842         float gedep = NAN;
0843 
0844         float efromtruth = NAN;
0845 
0846         if (primary)
0847         {
0848           gparticleID = primary->get_track_id();
0849           gflavor = primary->get_pid();
0850 
0851           PHG4Shower* shower = trutheval->get_primary_shower(primary);
0852           if (shower)
0853           {
0854             gnhits = shower->get_nhits(trutheval->get_caloid());
0855           }
0856           else
0857           {
0858             gnhits = 0.0;
0859           }
0860           gpx = primary->get_px();
0861           gpy = primary->get_py();
0862           gpz = primary->get_pz();
0863           ge = primary->get_e();
0864 
0865           gpt = std::sqrt(gpx * gpx + gpy * gpy);
0866           if (gpt != 0.0)
0867           {
0868             geta = std::asinh(gpz / gpt);
0869           }
0870           gphi = std::atan2(gpy, gpx);
0871 
0872           PHG4VtxPoint* vtx = trutheval->get_vertex(primary);
0873 
0874           if (vtx)
0875           {
0876             gvx = vtx->get_x();
0877             gvy = vtx->get_y();
0878             gvz = vtx->get_z();
0879           }
0880 
0881           gembed = trutheval->get_embed(primary);
0882           gedep = trutheval->get_shower_energy_deposit(primary);
0883 
0884           efromtruth = towereval->get_energy_contribution(tower, primary);
0885         }
0886 
0887         float tower_data[] = {(float) _ievent,
0888                               towerid,
0889                               ieta,
0890                               iphi,
0891                               eta,
0892                               phi,
0893                               e,
0894                               x,
0895                               y,
0896                               z,
0897                               gparticleID,
0898                               gflavor,
0899                               gnhits,
0900                               geta,
0901                               gphi,
0902                               ge,
0903                               gpt,
0904                               gvx,
0905                               gvy,
0906                               gvz,
0907                               gembed,
0908                               gedep,
0909                               efromtruth};
0910 
0911         _ntp_tower->Fill(tower_data);
0912         _tower_debug->Fill();  // Added by Barak (see above for explanation)
0913       }
0914     }
0915   }
0916   
0917 
0918   //------------------------
0919   // fill the Cluster NTuple
0920   //------------------------
0921 
0922   if (_do_cluster_eval)
0923   {
0924     if (Verbosity() > 1)
0925     {
0926       std::cout << "CaloEvaluator::filling gcluster ntuple..." << std::endl;
0927     }
0928 
0929     GlobalVertexMap* vertexmap = findNode::getClass<GlobalVertexMap>(topNode, "GlobalVertexMap");
0930 
0931     std::string clusternode = "CLUSTER_" + _caloname;
0932     if (_use_towerinfo)
0933     {
0934       clusternode = "CLUSTER_CALIB_" + _caloname;
0935     }
0936     RawClusterContainer* clusters = findNode::getClass<RawClusterContainer>(topNode, clusternode.c_str());
0937     if (!clusters)
0938     {
0939       std::cout << PHWHERE << " ERROR: Can't find " << clusternode << std::endl;
0940       exit(-1);
0941     }
0942 
0943     // for every cluster
0944 
0945     for (const auto& iterator : clusters->getClustersMap())
0946     {
0947       RawCluster* cluster = iterator.second;
0948 
0949       //    for (unsigned int icluster = 0; icluster < clusters->size(); icluster++)
0950       //    {
0951       //      RawCluster* cluster = clusters->getCluster(icluster);
0952 
0953       if (cluster->get_energy() < _reco_e_threshold)
0954       {
0955         continue;
0956       }
0957 
0958       float clusterID = cluster->get_id();
0959       float ntowers = cluster->getNTowers();
0960       float x = cluster->get_x();
0961       float y = cluster->get_y();
0962       float z = cluster->get_z();
0963       float eta = NAN;
0964       float phi = cluster->get_phi();
0965       float e = cluster->get_energy();
0966 
0967       // require vertex for cluster eta calculation
0968       if (vertexmap)
0969       {
0970         if (!vertexmap->empty())
0971         {
0972           GlobalVertex* vertex = (vertexmap->begin()->second);
0973 
0974           eta =
0975               RawClusterUtility::GetPseudorapidity(
0976                   *cluster,
0977                   CLHEP::Hep3Vector(vertex->get_x(), vertex->get_y(), vertex->get_z()));
0978         }
0979       }
0980 
0981       PHG4Particle* primary = clustereval->max_truth_primary_particle_by_energy(cluster);
0982 
0983       float gparticleID = NAN;
0984       float gflavor = NAN;
0985 
0986       float gnhits = NAN;
0987       float gpx = NAN;
0988       float gpy = NAN;
0989       float gpz = NAN;
0990       float ge = NAN;
0991 
0992       float gpt = NAN;
0993       float geta = NAN;
0994       float gphi = NAN;
0995 
0996       float gvx = NAN;
0997       float gvy = NAN;
0998       float gvz = NAN;
0999 
1000       float gembed = NAN;
1001       float gedep = NAN;
1002 
1003       float efromtruth = NAN;
1004 
1005       if (primary)
1006       {
1007         gparticleID = primary->get_track_id();
1008         gflavor = primary->get_pid();
1009 
1010         PHG4Shower* shower = trutheval->get_primary_shower(primary);
1011         if (shower)
1012         {
1013           gnhits = shower->get_nhits(trutheval->get_caloid());
1014         }
1015         else
1016         {
1017           gnhits = 0.0;
1018         }
1019         gpx = primary->get_px();
1020         gpy = primary->get_py();
1021         gpz = primary->get_pz();
1022         ge = primary->get_e();
1023 
1024         gpt = std::sqrt(gpx * gpx + gpy * gpy);
1025         if (gpt != 0.0)
1026         {
1027           geta = std::asinh(gpz / gpt);
1028         }
1029         gphi = std::atan2(gpy, gpx);
1030 
1031         PHG4VtxPoint* vtx = trutheval->get_vertex(primary);
1032 
1033         if (vtx)
1034         {
1035           gvx = vtx->get_x();
1036           gvy = vtx->get_y();
1037           gvz = vtx->get_z();
1038         }
1039 
1040         gembed = trutheval->get_embed(primary);
1041         gedep = trutheval->get_shower_energy_deposit(primary);
1042 
1043         efromtruth = clustereval->get_energy_contribution(cluster,
1044                                                           primary);
1045       }
1046 
1047       float cluster_data[] = {(float) _ievent,
1048                               clusterID,
1049                               ntowers,
1050                               eta,
1051                               x,
1052                               y,
1053                               z,
1054                               phi,
1055                               e,
1056                               gparticleID,
1057                               gflavor,
1058                               gnhits,
1059                               geta,
1060                               gphi,
1061                               ge,
1062                               gpt,
1063                               gvx,
1064                               gvy,
1065                               gvz,
1066                               gembed,
1067                               gedep,
1068                               efromtruth};
1069 
1070       _ntp_cluster->Fill(cluster_data);
1071     }
1072   }
1073 
1074   return;
1075 }