File indexing completed on 2025-12-16 09:20:27
0001
0002
0003
0004
0005
0006
0007 #include "TpcDirectLaserReconstruction.h"
0008
0009 #include "TpcSpaceChargeMatrixContainerv1.h"
0010
0011 #include <g4detectors/PHG4TpcGeom.h>
0012 #include <g4detectors/PHG4TpcGeomContainer.h>
0013
0014 #include <trackbase/ActsTrackingGeometry.h>
0015 #include <trackbase/TpcDefs.h>
0016 #include <trackbase/TrkrHitSet.h>
0017 #include <trackbase/TrkrHitSetContainer.h>
0018 #include <trackbase/TrkrHitv2.h> // for TrkrHit
0019 #include <trackbase_historic/SvtxTrack.h>
0020 #include <trackbase_historic/SvtxTrackMap.h>
0021 #include <trackbase_historic/SvtxTrackState_v1.h>
0022
0023 #include <fun4all/Fun4AllReturnCodes.h>
0024
0025 #include <phool/getClass.h>
0026
0027 #include <TFile.h>
0028 #include <TH1.h>
0029 #include <TH2.h>
0030 #include <TH3.h>
0031 #include <TNtuple.h>
0032 #include <TVector3.h>
0033
0034 #include <algorithm>
0035 #include <boost/format.hpp>
0036
0037 #include <cassert>
0038 #include <cmath>
0039
0040 namespace
0041 {
0042
0043
0044 template <class T>
0045 class range_adaptor
0046 {
0047 public:
0048 explicit range_adaptor(const T& range)
0049 : m_range(range)
0050 {
0051 }
0052 inline const typename T::first_type& begin() { return m_range.first; }
0053 inline const typename T::second_type& end() { return m_range.second; }
0054
0055 private:
0056 T m_range;
0057 };
0058
0059
0060 template <class T>
0061 inline constexpr T square(const T& x)
0062 {
0063 return x * x;
0064 }
0065
0066
0067 template <class T>
0068 inline constexpr T get_r(const T& x, const T& y)
0069 {
0070 return std::sqrt(square(x) + square(y));
0071 }
0072
0073
0074 std::pair<double, double> line_circle_intersection(const TVector3& p, const TVector3& d, double radius)
0075 {
0076 std::pair<double, double> ret = std::make_pair(-1, -1);
0077
0078 const double A = square(d.x()) + square(d.y());
0079 const double B = 2 * p.x() * d.x() + 2 * p.y() * d.y();
0080 const double C = square(p.x()) + square(p.y()) - square(radius);
0081 const double delta = square(B) - 4 * A * C;
0082 if (delta < 0)
0083 {
0084 return ret;
0085 }
0086
0087
0088 const double tup = (-B + std::sqrt(delta)) / (2 * A);
0089 const double tdn = (-B - sqrt(delta)) / (2 * A);
0090
0091 ret = std::make_pair(tup, tdn);
0092
0093 return ret;
0094
0095
0096
0097
0098
0099
0100
0101
0102
0103
0104 }
0105
0106
0107 inline std::ostream& operator<<(std::ostream& out, const TVector3& vector)
0108 {
0109 out << "( " << vector.x() << ", " << vector.y() << ", " << vector.z() << ")";
0110 return out;
0111 }
0112
0113
0114 template <class T>
0115 inline constexpr T delta_phi(const T& phi)
0116 {
0117 if (phi >= M_PI)
0118 {
0119 return phi - 2 * M_PI;
0120 }
0121 else if (phi < -M_PI)
0122 {
0123 return phi + 2 * M_PI;
0124 }
0125 else
0126 {
0127 return phi;
0128 }
0129 }
0130
0131
0132 static constexpr float m_phimin = 0;
0133 static constexpr float m_phimax = 2. * M_PI;
0134
0135
0136
0137 static constexpr float m_rmin = 20;
0138 static constexpr float m_rmax = 78;
0139
0140 }
0141
0142
0143 TpcDirectLaserReconstruction::TpcDirectLaserReconstruction(const std::string& name)
0144 : SubsysReco(name)
0145 , PHParameterInterface(name)
0146 , m_matrix_container(new TpcSpaceChargeMatrixContainerv1)
0147 {
0148 InitializeParameters();
0149 }
0150
0151
0152 int TpcDirectLaserReconstruction::Init(PHCompositeNode* )
0153 {
0154 m_total_hits = 0;
0155 m_matched_hits = 0;
0156 m_accepted_clusters = 0;
0157
0158 if (m_savehistograms)
0159 {
0160 create_histograms();
0161 }
0162 return Fun4AllReturnCodes::EVENT_OK;
0163 }
0164
0165
0166 int TpcDirectLaserReconstruction::InitRun(PHCompositeNode* )
0167 {
0168 UpdateParametersWithMacro();
0169 m_max_dca = get_double_param("directlaser_max_dca");
0170 m_max_drphi = get_double_param("directlaser_max_drphi");
0171 m_max_dz = get_double_param("directlaser_max_dz");
0172
0173
0174 if (Verbosity())
0175 {
0176 std::cout
0177 << "TpcDirectLaserReconstruction::InitRun\n"
0178 << " m_outputfile: " << m_outputfile << "\n"
0179 << " m_max_dca: " << m_max_dca << "\n"
0180 << " m_max_drphi: " << m_max_drphi << "\n"
0181 << " m_max_dz: " << m_max_dz << "\n"
0182 << std::endl;
0183
0184
0185 m_matrix_container->identify();
0186 }
0187
0188 return Fun4AllReturnCodes::EVENT_OK;
0189 }
0190
0191
0192 int TpcDirectLaserReconstruction::process_event(PHCompositeNode* topNode)
0193 {
0194
0195 const auto res = load_nodes(topNode);
0196 if (res != Fun4AllReturnCodes::EVENT_OK)
0197 {
0198 return res;
0199 }
0200
0201 m_zmax = m_tGeometry->get_max_driftlength() + m_tGeometry->get_CM_halfwidth();
0202 m_zmin = -m_zmax;
0203
0204 process_tracks();
0205 return Fun4AllReturnCodes::EVENT_OK;
0206 }
0207
0208
0209 int TpcDirectLaserReconstruction::End(PHCompositeNode* )
0210 {
0211
0212 if (m_matrix_container)
0213 {
0214 std::unique_ptr<TFile> outputfile(TFile::Open(m_outputfile.c_str(), "RECREATE"));
0215 outputfile->cd();
0216 m_matrix_container->Write("TpcSpaceChargeMatrixContainer");
0217 }
0218
0219
0220 if (m_savehistograms && m_histogramfile)
0221 {
0222 m_histogramfile->cd();
0223
0224
0225 for (const auto& o : std::initializer_list<TObject*>({h_dca_layer, h_deltarphi_layer_south, h_deltarphi_layer_north, h_deltaz_layer, h_deltar_r, h_deltheta_delphi, h_deltheta_delphi_1, h_deltheta_delphi_2, h_deltheta_delphi_3, h_deltheta_delphi_4, h_deltheta_delphi_5, h_deltheta_delphi_6, h_deltheta_delphi_7, h_deltheta_delphi_8, h_entries, h_hits, h_hits_reco, h_adc, h_adc_reco, h_adc_sum, h_adc_sum_reco, h_adc_sum_ratio, h_adc_sum_ratio_true, h_adc_vs_DCA_true, h_adc_sum_ratio_lasrangle, h_num_sum, h_num_sum_reco, h_num_sum_ratio,
0226 h_num_sum_ratio_true, h_adc_sum_ratio_true, h_adc_sum_ratio, h_num_sum_ratio_lasrangle, h_bright_hits_laser1, h_bright_hits_laser2, h_bright_hits_laser3, h_bright_hits_laser4, h_relangle_lasrangle,
0227 h_relangle_theta_lasrangle, h_relangle_phi_lasrangle, h_GEMs_hit, h_layers_hit, h_xy, h_xz, h_xy_pca, h_xz_pca, h_dca_path, h_zr, h_zr_pca, h_dz_z}))
0228
0229 {
0230 if (o)
0231 {
0232 o->Write();
0233 }
0234 }
0235 m_histogramfile->Close();
0236 }
0237
0238
0239
0240
0241
0242
0243
0244
0245
0246
0247
0248
0249 std::cout << "TpcDirectLaserReconstruction::End - m_total_hits: " << m_total_hits << std::endl;
0250 std::cout << "TpcDirectLaserReconstruction::End - m_matched_hits: " << m_matched_hits << std::endl;
0251 std::cout << "TpcDirectLaserReconstruction::End - m_accepted_clusters: " << m_accepted_clusters << std::endl;
0252 std::cout << "TpcDirectLaserReconstruction::End - fraction cluster/hits: " << m_accepted_clusters / m_total_hits << std::endl;
0253
0254 return Fun4AllReturnCodes::EVENT_OK;
0255 }
0256
0257
0258 void TpcDirectLaserReconstruction::SetDefaultParameters()
0259 {
0260
0261 set_default_double_param("directlaser_max_dca", 1.0);
0262
0263
0264
0265
0266
0267 set_default_double_param("directlaser_max_drphi", 2.);
0268 set_default_double_param("directlaser_max_dz", 2.);
0269 }
0270
0271
0272 void TpcDirectLaserReconstruction::set_grid_dimensions(int phibins, int rbins, int zbins)
0273 {
0274 m_matrix_container->set_grid_dimensions(phibins, rbins, zbins);
0275 }
0276
0277
0278 int TpcDirectLaserReconstruction::load_nodes(PHCompositeNode* topNode)
0279 {
0280 m_geom_container = findNode::getClass<PHG4TpcGeomContainer>(topNode, "TPCGEOMCONTAINER");
0281 assert(m_geom_container);
0282
0283
0284 m_tGeometry = findNode::getClass<ActsGeometry>(topNode, "ActsGeometry");
0285 assert(m_tGeometry);
0286
0287
0288 m_track_map = findNode::getClass<SvtxTrackMap>(topNode, "SvtxTrackMap");
0289 assert(m_track_map);
0290
0291
0292 m_hit_map = findNode::getClass<TrkrHitSetContainer>(topNode, "TRKR_HITSET");
0293 assert(m_hit_map);
0294
0295 return Fun4AllReturnCodes::EVENT_OK;
0296 }
0297
0298
0299 void TpcDirectLaserReconstruction::create_histograms()
0300 {
0301 std::cout << "TpcDirectLaserReconstruction::makeHistograms - writing evaluation histograms to: " << m_histogramfilename << std::endl;
0302 m_histogramfile.reset(new TFile(m_histogramfilename.c_str(), "RECREATE"));
0303 m_histogramfile->cd();
0304
0305
0306 h_dca_layer = new TH2F("dca_layer", ";radius; DCA (cm)", 78, 0, 78, 500, 0, 20);
0307 h_deltarphi_layer_north = new TH2F("deltarphi_layer_north", ";radius; r.#Delta#phi_{track-cluster} (cm)", 78, 0, 78, 2000, -5, 5);
0308 h_deltarphi_layer_south = new TH2F("deltarphi_layer_south", ";radius; r.#Delta#phi_{track-cluster} (cm)", 78, 0, 78, 2000, -5, 5);
0309 h_deltaz_layer = new TH2F("deltaz_layer", ";radius; #Deltaz_{track-cluster} (cm)", 78, 0, 78, 2000, -20, 20);
0310 h_deltar_r = new TH2F("deltar_r", ";radius;#Deltar_{track-cluster} (cm)", 78, 0, 78, 2000, -3, 3);
0311
0312 h_xy = new TH2F("h_xy", " x vs y", 320, -80, 80, 320, -80, 80);
0313 h_xz = new TH2F("h_xz", " x vs z", 320, -80, 80, 440, -110, 110);
0314 h_xy_pca = new TH2F("h_xy_pca", " x vs y pca", 320, -80, 80, 320, -80, 80);
0315 h_xz_pca = new TH2F("h_xz_pca", " x vs z pca", 320, -80, 80, 440, -110, 110);
0316 h_dca_path = new TH2F("h_dca_path", " dca vs pathlength", 440, 0, 110, 100, 0, 20);
0317 h_zr = new TH2F("h_zr", " z vs r", 440, -110, 110, 1000, 28, 80);
0318 h_zr->GetXaxis()->SetTitle("z");
0319 h_zr->GetYaxis()->SetTitle("rad");
0320 h_zr_pca = new TH2F("h_zr_pca", " z vs r pca", 440, -110, 110, 1000, 28, 80);
0321 h_dz_z = new TH2F("h_dz_z", " dz vs z", 440, -110, 110, 1000, -20, 20);
0322
0323 h_hits = new TNtuple("hits", "raw hits", "x:y:z:adc");
0324 h_hits_reco = new TNtuple("hits_reco", "raw hits", "x:y:z:adc");
0325
0326 h_adc = new TH1F("adc", "pedestal subtracted ADC spectra of ALL lasers (MCTRUTH direction)", 1063, -19.5, 1042.5);
0327 h_adc->SetXTitle("ADC - PEDESTAL [ADU]");
0328
0329 h_adc_reco = new TH1F("adc_reco", "pedestal subtracted ADC spectra of ALL lasers (RECO DIRECTION)", 1063, -19.5, 1042.5);
0330 h_adc_reco->SetXTitle("ADC - PEDESTAL [ADU]");
0331
0332 h_adc_sum = new TH1F("adc_sum", "non-pedestal subtracted sum of ADC for each laser (MCTRUTH DIRECTION)", 1600, -0.5, 159999.5);
0333 h_adc_sum->SetXTitle("#Sigma ADC [ADU]");
0334
0335 h_adc_sum_reco = new TH1F("adc_sum_reco", "non-pedestal subtracted sum of ADC for each laser (RECO DIRECTION)", 1600, -0.5, 159999.5);
0336 h_adc_sum_reco->SetXTitle("#Sigma ADC [ADU]");
0337
0338 h_num_sum = new TH1F("num_sum", "sum of hits for each laser (MCTRUTH DIRECTION)", 1600, -0.5, 7999.5);
0339 h_num_sum->SetXTitle("#Sigma N_{hits}^{laser X} [arb. units]");
0340
0341 h_num_sum_reco = new TH1F("num_sum_reco", "sum of hits for each laser (RECO DIRECTION)", 1600, -0.5, 7999.5);
0342 h_num_sum_reco->SetXTitle("#Sigma N_{hits}^{laser X} [arb. units]");
0343
0344 h_adc_sum_ratio = new TH1F("adc_sum_ratio", "RATIO of non-pedestal subtracted sum of ADC for each laser (RECO DIRECTION/MCTRUTH DIRECTION)", 120, -0.5, 5.5);
0345 h_adc_sum_ratio->SetXTitle("#frac{#Sigma^{RECO} ADC}{#Sigma^{MCTRUTH} ADC} [arb. units]");
0346
0347 h_adc_sum_ratio_true = new TH1F("adc_sum_ratio_true", "RATIO of non-pedestal subtracted sum of ADC for each laser (MCTRUTH DIRECTION/ALL )", 120, -0.5, 5.5);
0348 h_adc_sum_ratio_true->SetXTitle("#frac{#Sigma^{MCTRUTH} ADC}{#Sigma^{ALL} ADC} [arb. units]");
0349
0350 h_num_sum_ratio = new TH1F("num_sum_ratio", "RATIO of number of hits associated to each laser (RECO DIRECTION/MCTRUTH DIRECTION)", 120, -0.5, 5.5);
0351 h_num_sum_ratio->SetXTitle("#frac{#Sigma^{RECO} N_{hits}^{laser X}}{#Sigma^{MCTRUTH} N_{hits}^{laser X}} [arb. units]");
0352
0353 h_num_sum_ratio_true = new TH1F("num_sum_ratio_true", "RATIO of number of hits associated to each laser (MCTRUTH DIRECTION/ALL )", 120, -0.5, 5.5);
0354 h_num_sum_ratio_true->SetXTitle("#frac{#Sigma^{MCTRUTH} N_{hits}^{laser X}}{#Sigma^{ALL} N_{hits}^{laser X}} [arb. units]");
0355
0356 h_adc_vs_DCA_true = new TH2F("adc_vs_DCA_true", "PEDESTAL SUBTRACTED ADC vs DCA for ALL HITS from ALL LASERS", 100, 0, 100, 1024, -0.5, 1023.5);
0357 h_adc_vs_DCA_true->SetXTitle("DCA [mm]");
0358 h_adc_vs_DCA_true->SetYTitle("ADC - PEDESTAL [ADU]");
0359
0360 h_adc_sum_ratio_lasrangle = new TH2F("adc_sum_ratio_lasrangle", "RATIO of non-pedestal subtracted sum of ADC for LASER 1 ONLY (RECO DIRECTION/MCTRUTH DIRECTION)", 91, -0.5, 90.5, 361, -0.5, 360.5);
0361 h_adc_sum_ratio_lasrangle->SetXTitle("#theta_{laser} (deg.)");
0362 h_adc_sum_ratio_lasrangle->SetYTitle("#phi_{laser} (deg.)");
0363
0364 h_num_sum_ratio_lasrangle = new TH2F("num_sum_ratio_lasrangle", "RATIO of sum of hits for LASER 1 ONLY (RECO DIRECTION/MCTRUTH DIRECTION)", 91, -0.5, 90.5, 361, -0.5, 360.5);
0365 h_num_sum_ratio_lasrangle->SetXTitle("#theta_{laser} (deg.)");
0366 h_num_sum_ratio_lasrangle->SetYTitle("#phi_{laser} (deg.)");
0367
0368 h_bright_hits_laser1 = new TNtuple("bright_hits_laser1", "bright hits relative to laser 1", "x:y:z:deltheta:delphi");
0369 h_bright_hits_laser2 = new TNtuple("bright_hits_laser2", "bright hits relative to laser 2", "x:y:z:deltheta:delphi");
0370 h_bright_hits_laser3 = new TNtuple("bright_hits_laser3", "bright hits relative to laser 3", "x:y:z:deltheta:delphi");
0371 h_bright_hits_laser4 = new TNtuple("bright_hits_laser4", "bright hits relative to laser 4", "x:y:z:deltheta:delphi");
0372
0373
0374
0375
0376
0377
0378 h_relangle_lasrangle = new TH2F("relangle_lasrangle", "Difference in pointing angle from laser and relative angle histogram for ALL lasers", 102, -25.5, 25.5, 102, -25.5, 25.5);
0379 h_relangle_lasrangle->SetXTitle("#delta#theta (deg.)");
0380 h_relangle_lasrangle->SetYTitle("#delta#phi (deg.)");
0381
0382 h_relangle_theta_lasrangle = new TH2F("relangle_theta_lasrangle", "#delta#theta from laser and relative angle histogram for LASER 1 ONLY, ", 91, -0.5, 90.5, 361, -0.5, 360.5);
0383 h_relangle_theta_lasrangle->SetXTitle("#theta_{laser} (deg.)");
0384 h_relangle_theta_lasrangle->SetYTitle("#phi_{laser} (deg.)");
0385
0386 h_relangle_phi_lasrangle = new TH2F("relangle_phi_lasrangle", "#delta#phi from laser and relative angle histogram for LASER 1 ONLY", 91, -0.5, 90.5, 361, -0.5, 360.5);
0387 h_relangle_phi_lasrangle->SetXTitle("#theta_{laser} (deg.)");
0388 h_relangle_phi_lasrangle->SetYTitle("#phi_{laser} (deg.)");
0389
0390 h_deltheta_delphi = new TH2F("deltheta_delphi", "#Delta#theta, #Delta#phi for separation b/w TPC volume hits and ALL laser start points", 181, -0.5, 180.5, 361, -0.5, 360.5);
0391 h_deltheta_delphi->SetXTitle("#Delta#theta");
0392 h_deltheta_delphi->SetYTitle("#Delta#phi");
0393
0394 h_deltheta_delphi_1 = new TH2F("deltheta_delphi_1", "#Delta#theta, #Delta#phi for separation b/w TPC volume hits and LASER 1 only", 181, -0.5, 180.5, 361, -0.5, 360.5);
0395 h_deltheta_delphi_1->SetXTitle("#Delta#theta");
0396 h_deltheta_delphi_1->SetYTitle("#Delta#phi");
0397
0398 h_deltheta_delphi_2 = new TH2F("deltheta_delphi_2", "#Delta#theta, #Delta#phi for separation b/w TPC volume hits and LASER 2 only", 181, -0.5, 180.5, 361, -0.5, 360.5);
0399 h_deltheta_delphi_2->SetXTitle("#Delta#theta");
0400 h_deltheta_delphi_2->SetYTitle("#Delta#phi");
0401
0402 h_deltheta_delphi_3 = new TH2F("deltheta_delphi_3", "#Delta#theta, #Delta#phi for separation b/w TPC volume hits and LASER 3 only", 181, -0.5, 180.5, 361, -0.5, 360.5);
0403 h_deltheta_delphi_3->SetXTitle("#Delta#theta");
0404 h_deltheta_delphi_3->SetYTitle("#Delta#phi");
0405
0406 h_deltheta_delphi_4 = new TH2F("deltheta_delphi_4", "#Delta#theta, #Delta#phi for separation b/w TPC volume hits and LASER 4 only", 181, -0.5, 180.5, 361, -0.5, 360.5);
0407 h_deltheta_delphi_4->SetXTitle("#Delta#theta");
0408 h_deltheta_delphi_4->SetYTitle("#Delta#phi");
0409
0410 h_deltheta_delphi_5 = new TH2F("deltheta_delphi_5", "#Delta#theta, #Delta#phi for separation b/w TPC volume hits and LASER 5 only", 181, -0.5, 180.5, 361, -0.5, 360.5);
0411 h_deltheta_delphi_5->SetXTitle("#Delta#theta");
0412 h_deltheta_delphi_5->SetYTitle("#Delta#phi");
0413
0414 h_deltheta_delphi_6 = new TH2F("deltheta_delphi_6", "#Delta#theta, #Delta#phi for separation b/w TPC volume hits and LASER 6 only", 181, -0.5, 180.5, 361, -0.5, 360.5);
0415 h_deltheta_delphi_6->SetXTitle("#Delta#theta");
0416 h_deltheta_delphi_6->SetYTitle("#Delta#phi");
0417
0418 h_deltheta_delphi_7 = new TH2F("deltheta_delphi_7", "#Delta#theta, #Delta#phi for separation b/w TPC volume hits and LASER 7 only", 181, -0.5, 180.5, 361, -0.5, 360.5);
0419 h_deltheta_delphi_7->SetXTitle("#Delta#theta");
0420 h_deltheta_delphi_7->SetYTitle("#Delta#phi");
0421
0422 h_deltheta_delphi_8 = new TH2F("deltheta_delphi_8", "#Delta#theta, #Delta#phi for separation b/w TPC volume hits and LASER 8 only", 181, -0.5, 180.5, 361, -0.5, 360.5);
0423 h_deltheta_delphi_8->SetXTitle("#Delta#theta");
0424 h_deltheta_delphi_8->SetYTitle("#Delta#phi");
0425
0426 h_GEMs_hit = new TH1F("GEMS_hit", "Number of Unique GEM Modules hit for each laser", 8, 0, 8);
0427 h_layers_hit = new TH1F("layers_hit", "Number of Unique Layers hit for each laser", 8, 8, 8);
0428
0429 std::string GEM_bin_label;
0430 for (int GEMhistiter = 0; GEMhistiter < 8; GEMhistiter++)
0431 {
0432
0433 GEM_bin_label = (boost::format("laser %i") %(GEMhistiter + 1)).str();
0434 h_GEMs_hit->GetXaxis()->SetBinLabel(GEMhistiter + 1, GEM_bin_label.c_str());
0435 h_layers_hit->GetXaxis()->SetBinLabel(GEMhistiter + 1, GEM_bin_label.c_str());
0436 }
0437 h_GEMs_hit->SetYTitle("Number of Unique GEM Modules Hit");
0438 h_layers_hit->SetYTitle("Number of Unique Layers Hit");
0439
0440
0441
0442 if (m_matrix_container)
0443 {
0444 int phibins = 0;
0445 int rbins = 0;
0446 int zbins = 0;
0447 m_matrix_container->get_grid_dimensions(phibins, rbins, zbins);
0448 h_entries = new TH3F("entries", ";#phi;r (cm);z (cm)", phibins, m_phimin, m_phimax, rbins, m_rmin, m_rmax, zbins, m_zmin, m_zmax);
0449 }
0450 }
0451
0452
0453 void TpcDirectLaserReconstruction::process_tracks()
0454 {
0455 if (!(m_track_map && m_hit_map))
0456 {
0457 return;
0458 }
0459
0460
0461 for (auto& iter : *m_track_map)
0462 {
0463 process_track(iter.second);
0464 }
0465 }
0466
0467
0468 void TpcDirectLaserReconstruction::process_track(SvtxTrack* track)
0469 {
0470 std::multimap<unsigned int, std::pair<float, TVector3>> cluspos_map;
0471 std::set<unsigned int> layer_bin_set;
0472
0473
0474 const TVector3 origin(track->get_x(), track->get_y(), track->get_z());
0475 const TVector3 direction(track->get_px(), track->get_py(), track->get_pz());
0476 TVector3 direction2(track->get_px(), track->get_py(), track->get_pz());
0477
0478
0479
0480
0481
0482
0483
0484 const unsigned int trkid = track->get_id();
0485
0486 const float theta_orig = origin.Theta();
0487 const float phi_orig = origin.Phi();
0488
0489 float theta_trk = direction.Theta();
0490 direction2.RotateZ(-phi_orig);
0491 float phi_trk = direction2.Phi();
0492
0493 if (theta_trk > M_PI / 2.)
0494 {
0495 theta_trk = M_PI - theta_trk;
0496 }
0497
0498
0499
0500 while (phi_trk < m_phimin)
0501 {
0502 phi_trk += 2. * M_PI;
0503 }
0504 while (phi_trk >= m_phimax)
0505 {
0506 phi_trk -= 2. * M_PI;
0507 }
0508
0509
0510
0511
0512
0513
0514
0515 if (track->get_id() > 3)
0516 {
0517 phi_trk = (2 * M_PI) - phi_trk;
0518 }
0519
0520
0521
0522 std::cout << "---------- processing track " << track->get_id() << std::endl;
0523 std::cout << "TpcDirectLaserReconstruction::process_track - position: " << origin << " direction: " << direction << std::endl;
0524 std::cout << "Position Orientation - Theta: " << theta_orig * (180. / M_PI) << ", Phi: " << phi_orig * (180. / M_PI) << std::endl;
0525 std::cout << "Track Angle Direction - Theta: " << theta_trk * (180. / M_PI) << ", Phi: " << phi_trk * (180. / M_PI) << std::endl;
0526
0527
0528 int GEM_Mod_Arr[72] = {0};
0529
0530
0531 TrkrHitSetContainer::ConstRange hitsetrange = m_hit_map->getHitSets(TrkrDefs::TrkrId::tpcId);
0532
0533 float max_adc = 0.;
0534 float sum_adc_truth = 0;
0535 float sum_adc_truth_all = 0;
0536
0537 float sum_n_hits_truth = 0;
0538 float sum_n_hits_truth_all = 0;
0539
0540 for (auto hitsetitr = hitsetrange.first; hitsetitr != hitsetrange.second; ++hitsetitr)
0541 {
0542 const TrkrDefs::hitsetkey& hitsetkey = hitsetitr->first;
0543 const int side = TpcDefs::getSide(hitsetkey);
0544
0545 auto hitset = hitsetitr->second;
0546 const unsigned int layer = TrkrDefs::getLayer(hitsetkey);
0547 const auto layergeom = m_geom_container->GetLayerCellGeom(layer);
0548 const auto layer_center_radius = layergeom->get_radius();
0549
0550
0551
0552 const double AdcClockPeriod = layergeom->get_zstep();
0553 const unsigned short NTBins = (unsigned short) layergeom->get_zbins();
0554 const float tdriftmax = AdcClockPeriod * NTBins / 2.0;
0555
0556
0557 TrkrHitSet::ConstRange hitrangei = hitset->getHits();
0558
0559 for (auto hitr = hitrangei.first; hitr != hitrangei.second; ++hitr)
0560 {
0561 ++m_total_hits;
0562 sum_n_hits_truth_all++;
0563
0564 const unsigned short phibin = TpcDefs::getPad(hitr->first);
0565 const unsigned short zbin = TpcDefs::getTBin(hitr->first);
0566
0567 const double phi = layergeom->get_phicenter(phibin, side);
0568 const double x = layer_center_radius * cos(phi);
0569 const double y = layer_center_radius * sin(phi);
0570
0571 const double zdriftlength = layergeom->get_zcenter(zbin) * m_tGeometry->get_drift_velocity();
0572 double z = tdriftmax * m_tGeometry->get_drift_velocity() - zdriftlength;
0573 if (side == 0)
0574 {
0575 z *= -1;
0576 }
0577
0578 const TVector3 global(x, y, z);
0579
0580 float adc = (hitr->second->getAdc()) - m_pedestal;
0581 float adc_unsub = (hitr->second->getAdc());
0582 sum_adc_truth_all += adc_unsub;
0583
0584
0585
0586
0587
0588
0589 max_adc = std::max(adc, max_adc);
0590
0591
0592
0593 bool sameside = sameSign(global.z(), origin.z());
0594 const TVector3 oc(global.x() - origin.x(), global.y() - origin.y(), global.z() - origin.z());
0595 auto t = direction.Dot(oc) / square(direction.Mag());
0596 auto om = direction * t;
0597 const auto dca = (oc - om).Mag();
0598
0599 if (sameside)
0600 {
0601 h_adc_vs_DCA_true->Fill(dca, adc);
0602 }
0603
0604
0605
0606
0607
0608
0609
0610 TVector3 oc2(global.x() - origin.x(), global.y() - origin.y(), global.z() - origin.z());
0611
0612
0613
0614 oc2.RotateZ(-phi_orig);
0615
0616
0617
0618
0619
0620
0621
0622
0623
0624
0625 float deltheta = oc2.Theta();
0626 float delphi = oc2.Phi();
0627
0628 if (deltheta > M_PI / 2.)
0629 {
0630 deltheta = M_PI - deltheta;
0631 }
0632
0633 while (delphi < m_phimin)
0634 {
0635 delphi += 2. * M_PI;
0636 }
0637 while (delphi >= m_phimax)
0638 {
0639 delphi -= 2. * M_PI;
0640 }
0641
0642 if (track->get_id() > 3)
0643 {
0644 delphi = (2 * M_PI) - delphi;
0645 }
0646
0647
0648
0649 deltheta *= (180. / M_PI);
0650 delphi *= (180. / M_PI);
0651
0652
0653
0654
0655
0656
0657
0658
0659
0660
0661
0662
0663
0664
0665
0666
0667
0668
0669 if (sameside)
0670 {
0671
0672
0673 h_deltheta_delphi->Fill(deltheta, delphi);
0674
0675 if (track->get_id() == 0)
0676 {
0677
0678 h_deltheta_delphi_1->Fill(deltheta, delphi);
0679
0680
0681
0682
0683
0684
0685
0686
0687 }
0688 if (trkid == 1)
0689 {
0690 h_deltheta_delphi_2->Fill(deltheta, delphi);
0691
0692
0693
0694
0695
0696
0697 }
0698 if (trkid == 2)
0699 {
0700 h_deltheta_delphi_3->Fill(deltheta, delphi);
0701
0702
0703
0704
0705
0706
0707 }
0708 if (trkid == 3)
0709 {
0710 h_deltheta_delphi_4->Fill(deltheta, delphi);
0711
0712
0713
0714
0715
0716
0717 }
0718 if (trkid == 4)
0719 {
0720 h_deltheta_delphi_5->Fill(deltheta, delphi);
0721 }
0722 if (trkid == 5)
0723 {
0724 h_deltheta_delphi_6->Fill(deltheta, delphi);
0725 }
0726 if (trkid == 6)
0727 {
0728 h_deltheta_delphi_7->Fill(deltheta, delphi);
0729 }
0730 if (trkid == 7)
0731 {
0732 h_deltheta_delphi_8->Fill(deltheta, delphi);
0733 }
0734 }
0735
0736
0737 if (dca > m_max_dca)
0738 {
0739 continue;
0740 }
0741
0742 if (sameside)
0743 {
0744 sum_n_hits_truth++;
0745 h_adc->Fill(adc);
0746 sum_adc_truth += adc_unsub;
0747 }
0748
0749 ++m_matched_hits;
0750
0751
0752
0753
0754
0755
0756 float r2 = std::sqrt((x * x) + (y * y));
0757 float phi3 = phi;
0758 float z2 = z;
0759 while (phi3 < m_phimin)
0760 {
0761 phi3 += 2. * M_PI;
0762 }
0763 while (phi3 >= m_phimax)
0764 {
0765 phi3 -= 2. * M_PI;
0766 }
0767
0768 const int locateid = Locate(r2, phi3, z2);
0769
0770 if (z2 > m_zmin || z2 < m_zmax)
0771 {
0772 GEM_Mod_Arr[locateid - 1]++;
0773 }
0774
0775
0776 const auto cluspos_pair = std::make_pair(adc, global);
0777 cluspos_map.insert(std::make_pair(layer, cluspos_pair));
0778 layer_bin_set.insert(layer);
0779 }
0780 }
0781
0782 h_adc_sum->Fill(sum_adc_truth);
0783 h_num_sum->Fill(sum_n_hits_truth);
0784
0785 if (sum_adc_truth_all > 0)
0786 {
0787 h_adc_sum_ratio_true->Fill(sum_adc_truth / sum_adc_truth_all);
0788 }
0789
0790 if (sum_n_hits_truth_all > 0)
0791 {
0792 h_num_sum_ratio_true->Fill(sum_n_hits_truth / sum_n_hits_truth_all);
0793 }
0794
0795 int maxbin;
0796 int deltheta_max;
0797 int delphi_max;
0798 int dummy_z;
0799
0800 float theta_reco = 0;
0801 float phi_reco = 0;
0802
0803 float direction_reco;
0804 if (trkid < 4)
0805 {
0806 direction_reco = -1;
0807 }
0808 else
0809 {
0810 direction_reco = 1;
0811 }
0812
0813 TVector3 dir(0, 0, direction_reco);
0814
0815 if ((trkid == 0 && h_deltheta_delphi_1->GetEntries() > 0) && max_adc > 10)
0816 {
0817 h_deltheta_delphi_1->GetXaxis()->SetRangeUser(-5 + (theta_trk * (180. / M_PI)), 5 + (theta_trk * (180. / M_PI)));
0818 h_deltheta_delphi_1->GetYaxis()->SetRangeUser(-5 + (phi_trk * (180. / M_PI)), 5 + (phi_trk * (180. / M_PI)));
0819
0820 maxbin = h_deltheta_delphi_1->GetMaximumBin();
0821
0822 h_deltheta_delphi_1->GetBinXYZ(maxbin, deltheta_max, delphi_max, dummy_z);
0823
0824 std::cout << "Laser # " << (trkid + 1) << " naieve deltheta max: " << h_deltheta_delphi_1->GetXaxis()->GetBinCenter(deltheta_max) << ", naieve delphi max: " << h_deltheta_delphi_1->GetYaxis()->GetBinCenter(delphi_max) << std::endl;
0825
0826 h_relangle_lasrangle->Fill(h_deltheta_delphi_1->GetXaxis()->GetBinCenter(deltheta_max) - (theta_trk * (180. / M_PI)), h_deltheta_delphi_1->GetYaxis()->GetBinCenter(delphi_max) - (phi_trk * (180. / M_PI)));
0827
0828 h_relangle_theta_lasrangle->Fill(theta_trk * (180. / M_PI), phi_trk * (180. / M_PI), h_deltheta_delphi_1->GetXaxis()->GetBinCenter(deltheta_max) - (theta_trk * (180. / M_PI)));
0829 h_relangle_phi_lasrangle->Fill(theta_trk * (180. / M_PI), phi_trk * (180. / M_PI), h_deltheta_delphi_1->GetYaxis()->GetBinCenter(delphi_max) - (phi_trk * (180. / M_PI)));
0830
0831 theta_reco = h_deltheta_delphi_1->GetXaxis()->GetBinCenter(deltheta_max) * (M_PI / 180);
0832 phi_reco = h_deltheta_delphi_1->GetYaxis()->GetBinCenter(delphi_max) * (M_PI / 180);
0833 }
0834
0835 else if ((trkid == 1 && h_deltheta_delphi_2->GetEntries() > 0) && max_adc > 10)
0836 {
0837 h_deltheta_delphi_2->GetXaxis()->SetRangeUser(-5 + (theta_trk * (180. / M_PI)), 5 + (theta_trk * (180. / M_PI)));
0838 h_deltheta_delphi_2->GetYaxis()->SetRangeUser(-5 + (phi_trk * (180. / M_PI)), 5 + (phi_trk * (180. / M_PI)));
0839
0840 maxbin = h_deltheta_delphi_2->GetMaximumBin();
0841
0842 h_deltheta_delphi_2->GetBinXYZ(maxbin, deltheta_max, delphi_max, dummy_z);
0843
0844 std::cout << "Laser # " << (trkid + 1) << " naieve deltheta max: " << h_deltheta_delphi_2->GetXaxis()->GetBinCenter(deltheta_max) << ", naive delphi max: " << h_deltheta_delphi_2->GetYaxis()->GetBinCenter(delphi_max) << std::endl;
0845
0846 h_relangle_lasrangle->Fill(h_deltheta_delphi_2->GetXaxis()->GetBinCenter(deltheta_max) - (theta_trk * (180. / M_PI)), h_deltheta_delphi_2->GetYaxis()->GetBinCenter(delphi_max) - (phi_trk * (180. / M_PI)));
0847
0848
0849
0850
0851 theta_reco = h_deltheta_delphi_2->GetXaxis()->GetBinCenter(deltheta_max) * (M_PI / 180);
0852 phi_reco = h_deltheta_delphi_2->GetYaxis()->GetBinCenter(delphi_max) * (M_PI / 180);
0853 }
0854
0855 else if ((trkid == 2 && h_deltheta_delphi_3->GetEntries() > 0) && max_adc > 10)
0856 {
0857 h_deltheta_delphi_3->GetXaxis()->SetRangeUser(-5 + (theta_trk * (180. / M_PI)), 5 + (theta_trk * (180. / M_PI)));
0858 h_deltheta_delphi_3->GetYaxis()->SetRangeUser(-5 + (phi_trk * (180. / M_PI)), 5 + (phi_trk * (180. / M_PI)));
0859
0860 maxbin = h_deltheta_delphi_3->GetMaximumBin();
0861
0862 h_deltheta_delphi_3->GetBinXYZ(maxbin, deltheta_max, delphi_max, dummy_z);
0863
0864 std::cout << "Laser # " << (trkid + 1) << " naieve deltheta max: " << h_deltheta_delphi_3->GetXaxis()->GetBinCenter(deltheta_max) << ", naive delphi max: " << h_deltheta_delphi_3->GetYaxis()->GetBinCenter(delphi_max) << std::endl;
0865
0866 h_relangle_lasrangle->Fill(h_deltheta_delphi_3->GetXaxis()->GetBinCenter(deltheta_max) - (theta_trk * (180. / M_PI)), h_deltheta_delphi_3->GetYaxis()->GetBinCenter(delphi_max) - (phi_trk * (180. / M_PI)));
0867
0868
0869
0870
0871 theta_reco = h_deltheta_delphi_3->GetXaxis()->GetBinCenter(deltheta_max) * (M_PI / 180);
0872 phi_reco = h_deltheta_delphi_3->GetYaxis()->GetBinCenter(delphi_max) * (M_PI / 180);
0873 }
0874
0875 else if ((trkid == 3 && h_deltheta_delphi_4->GetEntries() > 0) && max_adc > 10)
0876 {
0877 h_deltheta_delphi_4->GetXaxis()->SetRangeUser(-5 + (theta_trk * (180. / M_PI)), 5 + (theta_trk * (180. / M_PI)));
0878 h_deltheta_delphi_4->GetYaxis()->SetRangeUser(-5 + (phi_trk * (180. / M_PI)), 5 + (phi_trk * (180. / M_PI)));
0879
0880 maxbin = h_deltheta_delphi_4->GetMaximumBin();
0881
0882 h_deltheta_delphi_4->GetBinXYZ(maxbin, deltheta_max, delphi_max, dummy_z);
0883
0884 std::cout << "Laser # " << (trkid + 1) << " naieve deltheta max: " << h_deltheta_delphi_4->GetXaxis()->GetBinCenter(deltheta_max) << ", naieve delphi max: " << h_deltheta_delphi_4->GetYaxis()->GetBinCenter(delphi_max) << std::endl;
0885
0886 h_relangle_lasrangle->Fill(h_deltheta_delphi_4->GetXaxis()->GetBinCenter(deltheta_max) - (theta_trk * (180. / M_PI)), h_deltheta_delphi_4->GetYaxis()->GetBinCenter(delphi_max) - (phi_trk * (180. / M_PI)));
0887
0888
0889
0890
0891 theta_reco = h_deltheta_delphi_4->GetXaxis()->GetBinCenter(deltheta_max) * (M_PI / 180);
0892 phi_reco = h_deltheta_delphi_4->GetYaxis()->GetBinCenter(delphi_max) * (M_PI / 180);
0893 }
0894
0895 else if ((trkid == 4 && h_deltheta_delphi_5->GetEntries() > 0) && max_adc > 10)
0896 {
0897 h_deltheta_delphi_5->GetXaxis()->SetRangeUser(-5 + (theta_trk * (180. / M_PI)), 5 + (theta_trk * (180. / M_PI)));
0898 h_deltheta_delphi_5->GetYaxis()->SetRangeUser(-5 + (phi_trk * (180. / M_PI)), 5 + (phi_trk * (180. / M_PI)));
0899
0900 maxbin = h_deltheta_delphi_5->GetMaximumBin();
0901
0902 h_deltheta_delphi_5->GetBinXYZ(maxbin, deltheta_max, delphi_max, dummy_z);
0903
0904 std::cout << "Laser # " << (trkid + 1) << " deltheta max: " << h_deltheta_delphi_5->GetXaxis()->GetBinCenter(deltheta_max) << ", delphi max: " << h_deltheta_delphi_5->GetYaxis()->GetBinCenter(delphi_max) << std::endl;
0905
0906 h_relangle_lasrangle->Fill(h_deltheta_delphi_5->GetXaxis()->GetBinCenter(deltheta_max) - (theta_trk * (180. / M_PI)), h_deltheta_delphi_5->GetYaxis()->GetBinCenter(delphi_max) - (phi_trk * (180. / M_PI)));
0907
0908
0909
0910
0911 theta_reco = h_deltheta_delphi_5->GetXaxis()->GetBinCenter(deltheta_max) * (M_PI / 180);
0912 phi_reco = h_deltheta_delphi_5->GetYaxis()->GetBinCenter(delphi_max) * (M_PI / 180);
0913 }
0914
0915 else if ((trkid == 5 && h_deltheta_delphi_6->GetEntries() > 0) && max_adc > 10)
0916 {
0917 h_deltheta_delphi_6->GetXaxis()->SetRangeUser(-5 + (theta_trk * (180. / M_PI)), 5 + (theta_trk * (180. / M_PI)));
0918 h_deltheta_delphi_6->GetYaxis()->SetRangeUser(-5 + (phi_trk * (180. / M_PI)), 5 + (phi_trk * (180. / M_PI)));
0919
0920 maxbin = h_deltheta_delphi_6->GetMaximumBin();
0921
0922 h_deltheta_delphi_6->GetBinXYZ(maxbin, deltheta_max, delphi_max, dummy_z);
0923
0924 std::cout << "Laser # " << (trkid + 1) << " deltheta max: " << h_deltheta_delphi_6->GetXaxis()->GetBinCenter(deltheta_max) << ", delphi max: " << h_deltheta_delphi_6->GetYaxis()->GetBinCenter(delphi_max) << std::endl;
0925
0926 h_relangle_lasrangle->Fill(h_deltheta_delphi_6->GetXaxis()->GetBinCenter(deltheta_max) - (theta_trk * (180. / M_PI)), h_deltheta_delphi_6->GetYaxis()->GetBinCenter(delphi_max) - (phi_trk * (180. / M_PI)));
0927
0928
0929
0930
0931 theta_reco = h_deltheta_delphi_6->GetXaxis()->GetBinCenter(deltheta_max) * (M_PI / 180);
0932 phi_reco = h_deltheta_delphi_6->GetYaxis()->GetBinCenter(delphi_max) * (M_PI / 180);
0933 }
0934
0935 else if ((trkid == 6 && h_deltheta_delphi_7->GetEntries() > 0) && max_adc > 10)
0936 {
0937 h_deltheta_delphi_7->GetXaxis()->SetRangeUser(-5 + (theta_trk * (180. / M_PI)), 5 + (theta_trk * (180. / M_PI)));
0938 h_deltheta_delphi_7->GetYaxis()->SetRangeUser(-5 + (phi_trk * (180. / M_PI)), 5 + (phi_trk * (180. / M_PI)));
0939
0940 maxbin = h_deltheta_delphi_7->GetMaximumBin();
0941
0942 h_deltheta_delphi_7->GetBinXYZ(maxbin, deltheta_max, delphi_max, dummy_z);
0943
0944 std::cout << "Laser # " << (trkid + 1) << " deltheta max: " << h_deltheta_delphi_7->GetXaxis()->GetBinCenter(deltheta_max) << ", delphi max: " << h_deltheta_delphi_7->GetYaxis()->GetBinCenter(delphi_max) << std::endl;
0945
0946 h_relangle_lasrangle->Fill(h_deltheta_delphi_7->GetXaxis()->GetBinCenter(deltheta_max) - (theta_trk * (180. / M_PI)), h_deltheta_delphi_7->GetYaxis()->GetBinCenter(delphi_max) - (phi_trk * (180. / M_PI)));
0947
0948
0949
0950
0951 theta_reco = h_deltheta_delphi_7->GetXaxis()->GetBinCenter(deltheta_max) * (M_PI / 180);
0952 phi_reco = h_deltheta_delphi_7->GetYaxis()->GetBinCenter(delphi_max) * (M_PI / 180);
0953 }
0954
0955 else if ((trkid == 7 && h_deltheta_delphi_8->GetEntries() > 0) && max_adc > 10)
0956 {
0957 h_deltheta_delphi_8->GetXaxis()->SetRangeUser(-5 + (theta_trk * (180. / M_PI)), 5 + (theta_trk * (180. / M_PI)));
0958 h_deltheta_delphi_8->GetYaxis()->SetRangeUser(-5 + (phi_trk * (180. / M_PI)), 5 + (phi_trk * (180. / M_PI)));
0959
0960 maxbin = h_deltheta_delphi_8->GetMaximumBin();
0961
0962 h_deltheta_delphi_8->GetBinXYZ(maxbin, deltheta_max, delphi_max, dummy_z);
0963
0964 std::cout << "Laser # " << (trkid + 1) << " deltheta max: " << h_deltheta_delphi_8->GetXaxis()->GetBinCenter(deltheta_max) << ", delphi max: " << h_deltheta_delphi_8->GetYaxis()->GetBinCenter(delphi_max) << std::endl;
0965
0966 h_relangle_lasrangle->Fill(h_deltheta_delphi_8->GetXaxis()->GetBinCenter(deltheta_max) - (theta_trk * (180. / M_PI)), h_deltheta_delphi_8->GetYaxis()->GetBinCenter(delphi_max) - (phi_trk * (180. / M_PI)));
0967
0968
0969
0970
0971 theta_reco = h_deltheta_delphi_8->GetXaxis()->GetBinCenter(deltheta_max) * (M_PI / 180);
0972 phi_reco = h_deltheta_delphi_8->GetYaxis()->GetBinCenter(delphi_max) * (M_PI / 180);
0973 }
0974
0975
0976 dir.RotateY(theta_reco * direction_reco);
0977 if (direction_reco == -1)
0978 {
0979 dir.RotateZ(phi_reco);
0980 }
0981 else
0982 {
0983 dir.RotateZ(-phi_reco);
0984 }
0985 dir.RotateZ(phi_orig);
0986
0987
0988
0989 TVector3 direction2_reco(dir.Px(), dir.Py(), dir.Pz());
0990
0991 float theta_trk_reco = dir.Theta();
0992 direction2_reco.RotateZ(-phi_orig);
0993 float phi_trk_reco = direction2_reco.Phi();
0994
0995 if (theta_trk_reco > M_PI / 2.)
0996 {
0997 theta_trk_reco = M_PI - theta_trk_reco;
0998 }
0999
1000
1001
1002 while (phi_trk_reco < m_phimin)
1003 {
1004 phi_trk_reco += 2. * M_PI;
1005 }
1006 while (phi_trk_reco >= m_phimax)
1007 {
1008 phi_trk_reco -= 2. * M_PI;
1009 }
1010
1011
1012
1013
1014
1015
1016
1017 if (track->get_id() > 3)
1018 {
1019 phi_trk_reco = (2 * M_PI) - phi_trk_reco;
1020 }
1021
1022 std::cout << "RECONSTRUCTED Track Angle Direction - Theta: " << theta_trk_reco * (180. / M_PI) << ", Phi: " << phi_trk_reco * (180. / M_PI) << std::endl;
1023
1024
1025
1026
1027
1028 float sum_adc_reco = 0;
1029 float sum_n_hits_reco = 0;
1030
1031 for (auto hitsetitr = hitsetrange.first; hitsetitr != hitsetrange.second; ++hitsetitr)
1032 {
1033 const TrkrDefs::hitsetkey& hitsetkey_2 = hitsetitr->first;
1034 const int side_2 = TpcDefs::getSide(hitsetkey_2);
1035
1036 auto hitset_2 = hitsetitr->second;
1037 const unsigned int layer_2 = TrkrDefs::getLayer(hitsetkey_2);
1038 const auto layergeom_2 = m_geom_container->GetLayerCellGeom(layer_2);
1039 const auto layer_center_radius_2 = layergeom_2->get_radius();
1040
1041
1042
1043 static constexpr double AdcClockPeriod_2 = 53.0;
1044 const unsigned short NTBins_2 = (unsigned short) layergeom_2->get_zbins();
1045 const float tdriftmax_2 = AdcClockPeriod_2 * NTBins_2 / 2.0;
1046
1047
1048 TrkrHitSet::ConstRange hitrangei_2 = hitset_2->getHits();
1049
1050 for (auto hitr = hitrangei_2.first; hitr != hitrangei_2.second; ++hitr)
1051 {
1052 const unsigned short phibin_2 = TpcDefs::getPad(hitr->first);
1053 const unsigned short zbin_2 = TpcDefs::getTBin(hitr->first);
1054
1055 const double phi_2 = layergeom_2->get_phicenter(phibin_2, side_2);
1056 const double x_2 = layer_center_radius_2 * cos(phi_2);
1057 const double y_2 = layer_center_radius_2 * sin(phi_2);
1058
1059 const double zdriftlength_2 = layergeom_2->get_zcenter(zbin_2) * m_tGeometry->get_drift_velocity();
1060 double z_2 = tdriftmax_2 * m_tGeometry->get_drift_velocity() - zdriftlength_2;
1061 if (side_2 == 0)
1062 {
1063 z_2 *= -1;
1064 }
1065
1066 const TVector3 global_2(x_2, y_2, z_2);
1067
1068 bool sameside_reco = sameSign(global_2.z(), origin.z());
1069
1070 float adc_2 = (hitr->second->getAdc()) - m_pedestal;
1071 float adc_2_unsub = (hitr->second->getAdc());
1072
1073
1074
1075 const TVector3 oc_2(global_2.x() - origin.x(), global_2.y() - origin.y(), global_2.z() - origin.z());
1076 auto t_2 = dir.Dot(oc_2) / square(dir.Mag());
1077 auto om_2 = dir * t_2;
1078 const auto dca_2 = (oc_2 - om_2).Mag();
1079
1080 if (dca_2 > m_max_dca)
1081 {
1082 continue;
1083 }
1084 if (sameside_reco)
1085 {
1086 sum_n_hits_reco++;
1087 h_adc_reco->Fill(adc_2);
1088 sum_adc_reco += adc_2_unsub;
1089 }
1090
1091
1092
1093
1094
1095
1096 }
1097 }
1098
1099 h_adc_sum_reco->Fill(sum_adc_reco);
1100 h_num_sum_reco->Fill(sum_n_hits_reco);
1101
1102 if (sum_adc_truth > 0)
1103 {
1104 h_adc_sum_ratio->Fill(sum_adc_reco / sum_adc_truth);
1105 if (trkid == 0)
1106 {
1107 h_adc_sum_ratio_lasrangle->Fill(theta_trk * (180. / M_PI), phi_trk * (180. / M_PI), (sum_adc_reco / sum_adc_truth));
1108 }
1109 }
1110
1111 if (sum_n_hits_truth > 0)
1112 {
1113 h_num_sum_ratio->Fill(sum_n_hits_reco / sum_n_hits_truth);
1114 if (trkid == 0)
1115 {
1116 h_num_sum_ratio_lasrangle->Fill(theta_trk * (180. / M_PI), phi_trk * (180. / M_PI), (sum_n_hits_reco / sum_n_hits_truth));
1117 }
1118 }
1119
1120 for (int GEMS_iter : GEM_Mod_Arr)
1121 {
1122 if (GEMS_iter > 0)
1123 {
1124 h_GEMs_hit->Fill(trkid + 0.5);
1125 }
1126 }
1127
1128 h_deltheta_delphi_1->Reset("ICESM");
1129 h_deltheta_delphi_2->Reset("ICESM");
1130 h_deltheta_delphi_3->Reset("ICESM");
1131 h_deltheta_delphi_4->Reset("ICESM");
1132 h_deltheta_delphi_5->Reset("ICESM");
1133 h_deltheta_delphi_6->Reset("ICESM");
1134 h_deltheta_delphi_7->Reset("ICESM");
1135 h_deltheta_delphi_8->Reset("ICESM");
1136
1137
1138
1139
1140
1141
1142
1143 for (auto layer : layer_bin_set)
1144 {
1145 PHG4TpcGeom* layergeom = m_geom_container->GetLayerCellGeom(layer);
1146 const auto layer_center_radius = layergeom->get_radius();
1147 const auto layer_inner_radius = layer_center_radius - layergeom->get_thickness() / 2.0;
1148 const auto layer_outer_radius = layer_center_radius + layergeom->get_thickness() / 2.0;
1149
1150 h_layers_hit->Fill(trkid + 0.5);
1151
1152
1153 auto tupdn = line_circle_intersection(origin, direction, layer_outer_radius);
1154 if (tupdn.first <= 0 && tupdn.second <= 0)
1155 {
1156 std::cout << " punt: layer " << layer << " layer outer radius " << layer_outer_radius << " tup " << tupdn.first << " tdn " << tupdn.second << std::endl;
1157 continue;
1158 }
1159 double layer_entry = tupdn.first;
1160 if (tupdn.second >= 0 && tupdn.second < tupdn.first)
1161 {
1162 layer_entry = tupdn.second;
1163 }
1164
1165 tupdn = line_circle_intersection(origin, direction, layer_inner_radius);
1166 if (tupdn.first <= 0 && tupdn.second <= 0)
1167 {
1168 std::cout << " punt: layer " << layer << " layer inner radius " << layer_inner_radius << " tup " << tupdn.first << " tdn " << tupdn.second << std::endl;
1169 continue;
1170 }
1171 double layer_exit = tupdn.first;
1172 if (tupdn.second > 0 && tupdn.second < tupdn.first)
1173 {
1174 layer_exit = tupdn.second;
1175 }
1176 if (Verbosity() > 2)
1177 {
1178 std::cout << " layer " << layer << " layer entry " << layer_entry << " layer exit " << layer_exit << std::endl;
1179 }
1180
1181
1182 tupdn = line_circle_intersection(origin, direction, layer_center_radius);
1183 double tup = tupdn.first;
1184 double tdn = tupdn.second;
1185
1186 double t = tup;
1187 if (tdn > 0 && tdn < tup)
1188 {
1189 t = tdn;
1190 }
1191 if (t < 0)
1192 {
1193 std::cout << " punt: layer " << layer << " layer center radius " << layer_center_radius << " t " << t << " tup " << tupdn.first << " tdn " << tupdn.second << std::endl;
1194 continue;
1195 }
1196
1197
1198 auto om = direction * t;
1199
1200
1201 const auto projection = origin + om;
1202
1203 double zmax = -999.;
1204 double zmin = 999.;
1205
1206 auto bin_residual = cluspos_map.equal_range(layer);
1207
1208 TVector3 clus_centroid(0, 0, 0);
1209 float wt = 0;
1210 for (auto mit = bin_residual.first; mit != bin_residual.second; ++mit)
1211 {
1212 float adc = (float) mit->second.first;
1213 TVector3 cluspos = mit->second.second;
1214
1215 if (fabs(cluspos.z() - projection.z()) > m_max_zrange)
1216 {
1217 continue;
1218 }
1219
1220 if (Verbosity() > 2)
1221 {
1222 std::cout << " layer " << layer << " adc " << adc << std::endl;
1223 std::cout << " cluspos " << cluspos.x() << " " << cluspos.y() << " " << cluspos.z()
1224 << " clus radius " << sqrt(square(cluspos.x()) + square(cluspos.y())) << std::endl;
1225 }
1226
1227 clus_centroid += cluspos * adc;
1228 wt += adc;
1229
1230 zmin = std::min(cluspos.z(), zmin);
1231 zmax = std::max(cluspos.z(), zmax);
1232 }
1233
1234 clus_centroid.SetX(clus_centroid.x() / wt);
1235 clus_centroid.SetY(clus_centroid.y() / wt);
1236 clus_centroid.SetZ(clus_centroid.z() / wt);
1237
1238 double zrange = zmax - zmin;
1239 if (zrange > m_max_zrange)
1240 {
1241 std::cout << " exceeded max zrange: zrange " << zrange << " max zrange " << m_max_zrange << std::endl;
1242 continue;
1243 }
1244
1245
1246
1247 const TVector3 oc(clus_centroid.x() - origin.x(), clus_centroid.y() - origin.y(), clus_centroid.z() - origin.z());
1248
1249 const auto dca = (oc - om).Mag();
1250
1251
1252 const auto pathlength = om.Mag();
1253
1254
1255 double ns_per_cm = 1e9 / 3e10;
1256 double dt = pathlength * ns_per_cm;
1257 double transit_dz = dt * m_tGeometry->get_drift_velocity();
1258 if (origin.z() > 0)
1259 {
1260 clus_centroid.SetZ(clus_centroid.z() + transit_dz);
1261 }
1262 else
1263 {
1264 clus_centroid.SetZ(clus_centroid.z() - transit_dz);
1265 }
1266
1267 if (Verbosity() > 0)
1268 {
1269 std::cout << " layer " << layer << " radius " << layer_center_radius << " wt " << wt << " t " << t << " dca " << dca
1270 << " pathlength " << pathlength << " transit_dz " << transit_dz << std::endl;
1271 std::cout << " clus_centroid " << clus_centroid.x() << " " << clus_centroid.y() << " " << clus_centroid.z() << " zrange " << zrange << std::endl;
1272 std::cout << " projection " << projection.x() << " " << projection.y() << " " << projection.z() << " dz " << clus_centroid.z() - projection.z() << std::endl;
1273 }
1274
1275
1276 SvtxTrackState_v1 state(pathlength);
1277 state.set_x(projection.x());
1278 state.set_y(projection.y());
1279 state.set_z(projection.z());
1280
1281 state.set_px(direction.x());
1282 state.set_py(direction.y());
1283 state.set_pz(direction.z());
1284 track->insert_state(&state);
1285
1286
1287
1288
1289
1290 const auto cluster_r = get_r(clus_centroid.x(), clus_centroid.y());
1291 const auto cluster_phi = std::atan2(clus_centroid.y(), clus_centroid.x());
1292 const auto cluster_z = clus_centroid.z();
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308 const auto cluster_rphi_error = 0.015;
1309 const auto cluster_z_error = 0.075;
1310
1311
1312 const auto track_phi = std::atan2(projection.y(), projection.x());
1313 const auto track_z = projection.z();
1314
1315
1316 const auto cosphi(std::cos(track_phi));
1317 const auto sinphi(std::sin(track_phi));
1318 const auto track_pphi = -state.get_px() * sinphi + state.get_py() * cosphi;
1319 const auto track_pr = state.get_px() * cosphi + state.get_py() * sinphi;
1320 const auto track_pz = state.get_pz();
1321 const auto talpha = -track_pphi / track_pr;
1322 const auto tbeta = -track_pz / track_pr;
1323
1324
1325 if (std::isnan(talpha))
1326 {
1327 std::cout << "TpcDirectLaserReconstruction::process_track - talpha is nan" << std::endl;
1328 continue;
1329 }
1330
1331 if (std::isnan(tbeta))
1332 {
1333 std::cout << "TpcDirectLaserReconstruction::process_track - tbeta is nan" << std::endl;
1334 continue;
1335 }
1336
1337
1338 const auto drp = cluster_r * delta_phi(cluster_phi - track_phi);
1339 const auto dz = cluster_z - track_z;
1340
1341
1342 if (std::isnan(drp))
1343 {
1344 std::cout << "TpcDirectLaserReconstruction::process_track - drp is nan" << std::endl;
1345 continue;
1346 }
1347
1348 if (std::isnan(dz))
1349 {
1350 std::cout << "TpcDirectLaserReconstruction::process_track - dz is nan" << std::endl;
1351 continue;
1352 }
1353
1354 if (m_savehistograms)
1355 {
1356 const float r = get_r(projection.x(), projection.y());
1357 const float dr = cluster_r - r;
1358 if (h_dca_layer)
1359 {
1360 h_dca_layer->Fill(r, dca);
1361 }
1362 if (h_deltarphi_layer_south && clus_centroid.z() < 0)
1363 {
1364 h_deltarphi_layer_south->Fill(r, drp);
1365 }
1366 if (h_deltarphi_layer_north && clus_centroid.z() > 0)
1367 {
1368 h_deltarphi_layer_north->Fill(r, drp);
1369 }
1370 if (h_deltaz_layer)
1371 {
1372 h_deltaz_layer->Fill(r, dz);
1373 }
1374 if (h_deltar_r)
1375 {
1376 h_deltar_r->Fill(r, dr);
1377 }
1378 if (h_entries)
1379 {
1380 auto phi = cluster_phi;
1381 while (phi < m_phimin)
1382 {
1383 phi += 2. * M_PI;
1384 }
1385 while (phi >= m_phimax)
1386 {
1387 phi -= 2. * M_PI;
1388 }
1389 h_entries->Fill(phi, cluster_r, cluster_z);
1390 }
1391 if (h_xy)
1392 {
1393 h_xy->Fill(clus_centroid.x(), clus_centroid.y());
1394 }
1395 if (h_xz)
1396 {
1397 h_xz->Fill(clus_centroid.x(), clus_centroid.z());
1398 }
1399 if (h_xy_pca)
1400 {
1401 h_xy_pca->Fill(projection.x(), projection.y());
1402 }
1403 if (h_xz_pca)
1404 {
1405 h_xz_pca->Fill(projection.x(), projection.z());
1406 }
1407 if (h_dca_path)
1408 {
1409 h_dca_path->Fill(pathlength, dca);
1410 }
1411 if (h_zr)
1412 {
1413 h_zr->Fill(clus_centroid.z(), cluster_r);
1414 }
1415 if (h_zr_pca)
1416 {
1417 h_zr_pca->Fill(projection.z(), r);
1418 }
1419 if (h_dz_z)
1420 {
1421 h_dz_z->Fill(projection.z(), clus_centroid.z() - projection.z());
1422 }
1423
1424 }
1425
1426
1427
1428
1429
1430
1431 const auto erp = square(cluster_rphi_error);
1432 const auto ez = square(cluster_z_error);
1433
1434
1435 if (std::isnan(erp))
1436 {
1437 std::cout << "TpcDirectLaserReconstruction::process_track - erp is nan" << std::endl;
1438 continue;
1439 }
1440
1441 if (std::isnan(ez))
1442 {
1443 std::cout << "TpcDirectLaserReconstruction::process_track - ez is nan" << std::endl;
1444 continue;
1445 }
1446
1447
1448 const auto i = get_cell_index(clus_centroid);
1449 if (i < 0)
1450 {
1451 if (Verbosity())
1452 {
1453 std::cout << "TpcDirectLaserReconstruction::process_track - invalid cell index"
1454 << " r: " << cluster_r
1455 << " phi: " << cluster_phi
1456 << " z: " << cluster_z
1457 << std::endl;
1458 }
1459 continue;
1460 }
1461
1462
1463
1464 m_matrix_container->add_to_lhs(i, 0, 0, 1. / erp);
1465 m_matrix_container->add_to_lhs(i, 0, 1, 0);
1466 m_matrix_container->add_to_lhs(i, 0, 2, talpha / erp);
1467
1468 m_matrix_container->add_to_lhs(i, 1, 0, 0);
1469 m_matrix_container->add_to_lhs(i, 1, 1, 1. / ez);
1470 m_matrix_container->add_to_lhs(i, 1, 2, tbeta / ez);
1471
1472 m_matrix_container->add_to_lhs(i, 2, 0, talpha / erp);
1473 m_matrix_container->add_to_lhs(i, 2, 1, tbeta / ez);
1474 m_matrix_container->add_to_lhs(i, 2, 2, square(talpha) / erp + square(tbeta) / ez);
1475
1476 m_matrix_container->add_to_rhs(i, 0, drp / erp);
1477 m_matrix_container->add_to_rhs(i, 1, dz / ez);
1478 m_matrix_container->add_to_rhs(i, 2, talpha * drp / erp + tbeta * dz / ez);
1479
1480
1481 m_matrix_container->add_to_entries(i);
1482
1483
1484 ++m_accepted_clusters;
1485 }
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496 }
1497
1498
1499 int TpcDirectLaserReconstruction::get_cell_index(const TVector3& global) const
1500 {
1501
1502 int phibins = 0;
1503 int rbins = 0;
1504 int zbins = 0;
1505 m_matrix_container->get_grid_dimensions(phibins, rbins, zbins);
1506
1507
1508
1509 float phi = std::atan2(global.y(), global.x());
1510 while (phi < m_phimin)
1511 {
1512 phi += 2. * M_PI;
1513 }
1514 while (phi >= m_phimax)
1515 {
1516 phi -= 2. * M_PI;
1517 }
1518 int iphi = phibins * (phi - m_phimin) / (m_phimax - m_phimin);
1519
1520
1521 const float r = get_r(global.x(), global.y());
1522 if (r < m_rmin || r >= m_rmax)
1523 {
1524 return -1;
1525 }
1526 int ir = rbins * (r - m_rmin) / (m_rmax - m_rmin);
1527
1528
1529 const float z = global.z();
1530 if (z < m_zmin || z >= m_zmax)
1531 {
1532 return -1;
1533 }
1534 int iz = zbins * (z - m_zmin) / (m_zmax - m_zmin);
1535
1536 return m_matrix_container->get_cell_index(iphi, ir, iz);
1537 }
1538
1539 int TpcDirectLaserReconstruction::Locate(float r, float phi, float z)
1540 {
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579 int GEM_Mod_ID;
1580
1581 const float Angle_Bins[13] = {23 * M_PI / 12, 1 * M_PI / 12, 3 * M_PI / 12, 5 * M_PI / 12, 7 * M_PI / 12, 9 * M_PI / 12, 11 * M_PI / 12, 13 * M_PI / 12, 15 * M_PI / 12, 17 * M_PI / 12, 19 * M_PI / 12, 21 * M_PI / 12, 23 * M_PI / 12};
1582
1583 const float r_bins[4] = {30, 46, 62, 78};
1584
1585 int angle_id = 0;
1586 int r_id = 0;
1587 int side_id = 0;
1588
1589 for (int r_iter = 0; r_iter < 3; r_iter++)
1590 {
1591 if (r > r_bins[r_iter] && r < r_bins[r_iter + 1])
1592 {
1593 r_id = r_iter + 1;
1594 break;
1595 }
1596 }
1597
1598 for (int ang_iter = 0; ang_iter < 12; ang_iter++)
1599 {
1600 if (phi > Angle_Bins[ang_iter] && phi < Angle_Bins[ang_iter + 1])
1601 {
1602 angle_id = ang_iter;
1603 break;
1604 }
1605 }
1606 if (z < 0)
1607 {
1608 side_id = 1;
1609 }
1610
1611 return GEM_Mod_ID = (36 * side_id) + (3 * angle_id) + r_id;
1612 }
1613
1614 float TpcDirectLaserReconstruction::GetRelPhi(float xorig, float yorig, float x, float y, float phiorig)
1615 {
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628 if (phiorig < 0)
1629 {
1630 phiorig += 2. * M_PI;
1631 }
1632 else if (phiorig > 2. * M_PI)
1633 {
1634 phiorig -= 2. * M_PI;
1635 }
1636
1637 float dx = x - xorig;
1638 float dy = y - yorig;
1639 float relphi = std::atan2(dy, dx) - phiorig;
1640 if (relphi < 0)
1641 {
1642 relphi += 2. * M_PI;
1643 }
1644 else if (relphi > 2. * M_PI)
1645 {
1646 relphi -= 2. * M_PI;
1647 }
1648
1649 return relphi;
1650 }
1651
1652 float TpcDirectLaserReconstruction::GetRelTheta(float xorig, float yorig, float zorig, float x, float y, float z, float thetaorig, float phiorig)
1653 {
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669 if (phiorig < 0)
1670 {
1671 phiorig += 2. * M_PI;
1672 }
1673 else if (phiorig > 2. * M_PI)
1674 {
1675 phiorig -= 2. * M_PI;
1676 }
1677
1678 if (thetaorig < 0)
1679 {
1680 thetaorig += 2. * M_PI;
1681 }
1682 else if (thetaorig > 2. * M_PI)
1683 {
1684 thetaorig -= 2. * M_PI;
1685 }
1686
1687 float dx = x - xorig;
1688 float dy = y - yorig;
1689 float dz = z - zorig;
1690 float r = std::sqrt(dx * dx + dy * dy + dz * dz);
1691
1692 float cos_beta = (dx * std::cos(phiorig) + dy * std::sin(phiorig)) / r;
1693 float sin_beta = dz / r;
1694
1695 float reltheta = std::acos(cos_beta * std::cos(thetaorig) + sin_beta * std::sin(thetaorig)) - M_PI / 2.;
1696 if (reltheta < 0)
1697 {
1698 reltheta += 2. * M_PI;
1699 }
1700 else if (reltheta > 2. * M_PI)
1701 {
1702 reltheta -= 2. * M_PI;
1703 }
1704
1705 return reltheta;
1706 }
1707
1708 bool TpcDirectLaserReconstruction::sameSign(float num1, float num2)
1709 {
1710 return (num1 >= 0 && num2 >= 0) || (num1 < 0 && num2 < 0);
1711 }