File indexing completed on 2025-08-06 08:18:04
0001
0002
0003
0004
0005
0006
0007
0008 #include "TpcCentralMembraneMatching.h"
0009
0010 #include <trackbase/CMFlashDifferenceContainerv1.h>
0011 #include <trackbase/CMFlashDifferencev1.h>
0012 #include <trackbase/LaserClusterContainerv1.h>
0013 #include <trackbase/LaserCluster.h>
0014 #include <trackbase/TpcDefs.h>
0015
0016 #include <ffaobjects/EventHeader.h>
0017
0018 #include <fun4all/Fun4AllReturnCodes.h>
0019
0020 #include <phool/PHCompositeNode.h>
0021 #include <phool/getClass.h>
0022 #include <phool/phool.h>
0023
0024 #include <TF1.h>
0025 #include <TFile.h>
0026 #include <TGraph.h>
0027 #include <TH1.h>
0028 #include <TH1F.h>
0029 #include <TH2.h>
0030 #include <TString.h>
0031 #include <TStyle.h>
0032 #include <TTree.h>
0033 #include <TVector3.h>
0034
0035 #include <boost/format.hpp>
0036
0037 #include <cmath>
0038 #include <iomanip>
0039 #include <set>
0040 #include <string>
0041
0042 int layerMins[3] = {16,23,39};
0043 int layerMaxes[3] = {22, 38, 54};
0044
0045 namespace
0046 {
0047 template <class T>
0048 inline constexpr T delta_phi(const T& phi)
0049 {
0050 if (phi > M_PI)
0051 {
0052 return phi - 2. * M_PI;
0053 }
0054 else if (phi <= -M_PI)
0055 {
0056 return phi + 2. * M_PI;
0057 }
0058 else
0059 {
0060 return phi;
0061 }
0062 }
0063
0064 template <class T>
0065 inline constexpr T square(const T& x)
0066 {
0067 return x * x;
0068 }
0069
0070 template <class T>
0071 inline T get_r(const T& x, const T& y)
0072 {
0073 return std::sqrt(square(x) + square(y));
0074 }
0075
0076
0077 [[maybe_unused]] std::ostream& operator<<(std::ostream& out, const Acts::Vector3& v)
0078 {
0079 out << "(" << v.x() << ", " << v.y() << ", " << v.z() << ")";
0080 return out;
0081 }
0082
0083
0084 [[maybe_unused]] void normalize_distortions(TpcDistortionCorrectionContainer* dcc)
0085 {
0086
0087 for (unsigned int i = 0; i < 2; ++i)
0088 {
0089
0090 for (int ip = 0; ip < dcc->m_hentries[i]->GetNbinsX(); ++ip)
0091 {
0092 for (int ir = 0; ir < dcc->m_hentries[i]->GetNbinsY(); ++ir)
0093 {
0094
0095 const auto entries = dcc->m_hentries[i]->GetBinContent(ip + 1, ir + 1);
0096 if (entries <= 1)
0097 {
0098 continue;
0099 }
0100
0101
0102 for (const auto& h : {dcc->m_hDRint[i], dcc->m_hDPint[i], dcc->m_hDZint[i]})
0103 {
0104 h->SetBinContent(ip + 1, ir + 1, h->GetBinContent(ip + 1, ir + 1) / entries);
0105 h->SetBinError(ip + 1, ir + 1, h->GetBinError(ip + 1, ir + 1) / entries);
0106 }
0107 }
0108 }
0109 }
0110 }
0111
0112
0113 [[maybe_unused]] void fill_guarding_bins(TpcDistortionCorrectionContainer* dcc)
0114 {
0115
0116 for (unsigned int i = 0; i < 2; ++i)
0117 {
0118 for (const auto& h : {dcc->m_hDRint[i], dcc->m_hDPint[i], dcc->m_hDZint[i], dcc->m_hentries[i]})
0119 {
0120
0121
0122
0123
0124
0125
0126 const auto phibins = h->GetNbinsX();
0127 const auto rbins = h->GetNbinsY();
0128 for (int ir = 0; ir < rbins; ++ir)
0129 {
0130
0131 h->SetBinContent(1, ir + 1, h->GetBinContent(phibins - 1, ir + 1));
0132 h->SetBinError(1, ir + 1, h->GetBinError(phibins - 1, ir + 1));
0133
0134
0135 h->SetBinContent(phibins, ir + 1, h->GetBinContent(2, ir + 1));
0136 h->SetBinError(phibins, ir + 1, h->GetBinError(2, ir + 1));
0137 }
0138
0139
0140 for (int iphi = 0; iphi < phibins; ++iphi)
0141 {
0142
0143 h->SetBinContent(iphi + 1, 1, h->GetBinContent(iphi + 1, 2));
0144 h->SetBinError(iphi + 1, 1, h->GetBinError(iphi + 1, 2));
0145
0146
0147 h->SetBinContent(iphi + 1, rbins, h->GetBinContent(iphi + 1, rbins - 1));
0148 h->SetBinError(iphi + 1, rbins, h->GetBinError(iphi + 1, rbins - 1));
0149 }
0150 }
0151 }
0152 }
0153
0154 }
0155
0156
0157 TpcCentralMembraneMatching::TpcCentralMembraneMatching(const std::string& name)
0158 : SubsysReco(name)
0159 {
0160
0161 CalculateCenters(nPads_R1, R1_e, nGoodStripes_R1_e, keepUntil_R1_e, nStripesIn_R1_e, nStripesBefore_R1_e, cx1_e, cy1_e);
0162 CalculateCenters(nPads_R1, R1, nGoodStripes_R1, keepUntil_R1, nStripesIn_R1, nStripesBefore_R1, cx1, cy1);
0163 CalculateCenters(nPads_R2, R2, nGoodStripes_R2, keepUntil_R2, nStripesIn_R2, nStripesBefore_R2, cx2, cy2);
0164 CalculateCenters(nPads_R3, R3, nGoodStripes_R3, keepUntil_R3, nStripesIn_R3, nStripesBefore_R3, cx3, cy3);
0165 }
0166
0167
0168 void TpcCentralMembraneMatching::set_grid_dimensions(int phibins, int rbins)
0169 {
0170 m_phibins = phibins;
0171 m_rbins = rbins;
0172 }
0173
0174
0175 double TpcCentralMembraneMatching::getPhiRotation_smoothed(TH1* hitHist, TH1* clustHist, bool side)
0176 {
0177
0178 hitHist->Smooth();
0179 clustHist->Smooth();
0180
0181
0182 TF1* f1 = new TF1(
0183 "f1", [&](double* x, double* p)
0184 { return p[0] * hitHist->GetBinContent(hitHist->FindBin((x[0] - p[1]) > M_PI ? x[0] - p[1] - 2 * M_PI : x[0] - p[1])); },
0185 -M_PI, M_PI, 2);
0186 f1->SetParNames("A", "shift");
0187 f1->SetParameters(1.0, 0.0);
0188 f1->SetParLimits(1, -M_PI / 18, M_PI / 18);
0189
0190 if (side && m_recoRotation[1][1] != -999)
0191 {
0192 f1->SetParameter(1, m_recoRotation[1][1]);
0193 }
0194 if (!side && m_recoRotation[0][1] != -999)
0195 {
0196 f1->SetParameter(1, m_recoRotation[0][1]);
0197 }
0198
0199 f1->SetNpx(500);
0200
0201 clustHist->Fit("f1", "IQ");
0202
0203
0204
0205
0206 return f1->GetParameter(1);
0207 }
0208
0209 void TpcCentralMembraneMatching::getRegionPhiRotation(bool side)
0210 {
0211 TH1* oddHist[3];
0212 TH1* evenHist[3];
0213
0214 TH1* oddTruthHist[3];
0215 TH1* evenTruthHist[3];
0216
0217 for (int i = 0; i < (int) m_reco_RPeaks[side].size(); i++)
0218 {
0219 int peak = reco_r_phi[side]->GetYaxis()->FindBin(m_reco_RPeaks[side][i]);
0220
0221 int region = 2;
0222 if (m_reco_RPeaks[side][i] < 41)
0223 {
0224 region = 0;
0225 }
0226 else if (m_reco_RPeaks[side][i] < 58)
0227 {
0228 region = 1;
0229 }
0230
0231 if (m_reco_RMatches[side][i] % 2)
0232 {
0233 if (!oddHist[region])
0234 {
0235 oddHist[region] = (TH1*) reco_r_phi[side]->ProjectionX((boost::format("oddR%d_%d") %region %side).str().c_str(), peak - 2, peak + 2);
0236 }
0237 else
0238 {
0239 oddHist[region]->Add(reco_r_phi[side]->ProjectionX((boost::format("oddR%d_%d") %region %side).str().c_str(), peak - 2, peak + 2));
0240 }
0241 }
0242 else
0243 {
0244 if (!evenHist[region])
0245 {
0246 evenHist[region] = (TH1*) reco_r_phi[side]->ProjectionX((boost::format("evenR%d_%d") %region %side).str().c_str(), peak - 2, peak + 2);
0247 }
0248 else
0249 {
0250 evenHist[region]->Add(reco_r_phi[side]->ProjectionX((boost::format("evenR%d_%d") %region %side).str().c_str(), peak - 2, peak + 2));
0251 }
0252 }
0253 }
0254
0255 for (int i = 0; i < (int) m_truth_RPeaks.size(); i++)
0256 {
0257 int peak = truth_r_phi[side]->GetYaxis()->FindBin(m_truth_RPeaks[i]);
0258
0259 int region = 2;
0260 if (m_truth_RPeaks[i] < 41)
0261 {
0262 region = 0;
0263 }
0264 else if (m_truth_RPeaks[i] < 58)
0265 {
0266 region = 1;
0267 }
0268
0269 if (i % 2)
0270 {
0271 if (!oddTruthHist[region])
0272 {
0273 oddTruthHist[region] = (TH1*) truth_r_phi[side]->ProjectionX((boost::format("oddTR%d_%d") %region %side).str().c_str(), peak - 1, peak + 1);
0274 }
0275 else
0276 {
0277 oddTruthHist[region]->Add(truth_r_phi[side]->ProjectionX((boost::format("oddTR%d_%d") %region %side).str().c_str(), peak - 1, peak + 1));
0278 }
0279 }
0280 else
0281 {
0282 if (!evenTruthHist[region])
0283 {
0284 evenTruthHist[region] = (TH1*) truth_r_phi[side]->ProjectionX((boost::format("evenTR%d_%d") %region %side).str().c_str(), peak - 1, peak + 1);
0285 }
0286 else
0287 {
0288 evenTruthHist[region]->Add(truth_r_phi[side]->ProjectionX((boost::format("evenTR%d_%d") %region %side).str().c_str(), peak - 1, peak + 1));
0289 }
0290 }
0291 }
0292
0293 int regions[3]{1, 0, 2};
0294
0295 double rotationEven[3]{-999, -999, -999};
0296 double rotationOdd[3]{-999, -999, -999};
0297
0298 for (int region : regions)
0299 {
0300 TF1* fEven = new TF1(
0301 "fEven", [&](double* x, double* p)
0302 { return p[0] * evenTruthHist[region]->GetBinContent(evenTruthHist[region]->FindBin((x[0] - p[1]) > M_PI ? x[0] - p[1] - 2 * M_PI : x[0] - p[1])); },
0303 -M_PI, M_PI, 2);
0304 fEven->SetParNames("A", "shift");
0305 fEven->SetParameters(1.0, 0.0);
0306 fEven->SetParLimits(1, -M_PI / 18, M_PI / 18);
0307 fEven->SetNpx(1000);
0308
0309 TF1* fOdd = new TF1(
0310 "fOdd", [&](double* x, double* p)
0311 { return p[0] * oddTruthHist[region]->GetBinContent(oddTruthHist[region]->FindBin((x[0] - p[1]) > M_PI ? x[0] - p[1] - 2 * M_PI : x[0] - p[1])); },
0312 -M_PI, M_PI, 2);
0313 fOdd->SetParNames("A", "shift");
0314 fOdd->SetParameters(1.0, 0.0);
0315 fOdd->SetParLimits(1, -M_PI / 18, M_PI / 18);
0316 fOdd->SetNpx(1000);
0317
0318 if (region != 1)
0319 {
0320 if (rotationEven[1] != -999)
0321 {
0322 fEven->SetParameter(1, rotationEven[1]);
0323 }
0324
0325 if (rotationOdd[1] != -999)
0326 {
0327 fOdd->SetParameter(1, rotationOdd[1]);
0328 }
0329 }
0330
0331 evenHist[region]->Fit("fEven", "IQ");
0332 oddHist[region]->Fit("fOdd", "IQ");
0333
0334 rotationEven[region] = fEven->GetParameter(1);
0335 rotationOdd[region] = fOdd->GetParameter(1);
0336
0337 m_recoRotation[side][region] = (rotationOdd[region] + rotationEven[region]) / 2.0;
0338 }
0339
0340
0341
0342
0343
0344
0345
0346
0347
0348
0349
0350
0351
0352
0353
0354
0355
0356
0357
0358
0359
0360
0361
0362
0363
0364
0365
0366
0367
0368
0369
0370
0371
0372
0373
0374
0375
0376
0377
0378
0379
0380
0381
0382
0383
0384
0385
0386
0387
0388
0389
0390
0391
0392
0393
0394
0395
0396
0397
0398
0399
0400
0401
0402
0403
0404
0405
0406
0407
0408
0409
0410
0411
0412
0413
0414
0415
0416
0417
0418
0419
0420
0421
0422
0423
0424
0425
0426
0427
0428
0429
0430
0431
0432
0433
0434
0435
0436
0437
0438
0439
0440
0441
0442
0443
0444
0445
0446
0447
0448
0449
0450
0451
0452
0453
0454
0455
0456
0457
0458
0459
0460
0461
0462
0463
0464
0465
0466
0467
0468
0469
0470
0471
0472
0473
0474
0475
0476
0477
0478
0479
0480
0481
0482
0483
0484
0485
0486
0487
0488
0489
0490
0491
0492
0493
0494
0495
0496
0497
0498
0499
0500
0501
0502
0503
0504
0505
0506
0507
0508
0509
0510
0511
0512
0513
0514
0515
0516
0517
0518
0519
0520
0521
0522
0523
0524
0525
0526
0527
0528
0529
0530
0531
0532
0533
0534
0535
0536
0537
0538
0539
0540
0541
0542
0543
0544
0545
0546
0547
0548
0549
0550
0551
0552
0553
0554
0555
0556
0557
0558
0559
0560
0561
0562
0563
0564
0565
0566
0567
0568
0569
0570
0571
0572
0573
0574
0575 }
0576
0577 std::vector<int> TpcCentralMembraneMatching::doGlobalRMatching(TH2* r_phi, bool side)
0578 {
0579 TH1D* proj = r_phi->ProjectionY("R_proj");
0580
0581 if (side)
0582 {
0583 m_m[1] = 0.0;
0584 m_b[1] = 0.0;
0585 }
0586 else
0587 {
0588 m_m[0] = 0.0;
0589 m_b[0] = 0.0;
0590 }
0591
0592 std::vector<double> rPeaks;
0593 std::vector<double> rHeights;
0594 std::vector<double> finalRPeaks;
0595 std::vector<double> finalRHeights;
0596 std::vector<std::vector<int>> groupR;
0597
0598 proj->GetXaxis()->SetRangeUser(0, 41);
0599 double maxR1 = proj->GetMaximum();
0600
0601 proj->GetXaxis()->SetRangeUser(41, 58);
0602 double maxR2 = proj->GetMaximum();
0603
0604 proj->GetXaxis()->SetRangeUser(58, 100);
0605 double maxR3 = proj->GetMaximum();
0606
0607 proj->GetXaxis()->SetRange(0, 0);
0608
0609 std::vector<int> hitMatches;
0610
0611 for (int i = 2; i < proj->GetNbinsX(); i++)
0612 {
0613
0614 if ((proj->GetBinCenter(i) < 41.0 && proj->GetBinContent(i) > 0.15 * maxR1) || (proj->GetBinCenter(i) >= 41.0 && proj->GetBinCenter(i) < 58.0 && proj->GetBinContent(i) > 0.15 * maxR2) || (proj->GetBinCenter(i) >= 58.0 && proj->GetBinContent(i) > 0.15 * maxR3))
0615 {
0616 rPeaks.push_back(proj->GetBinCenter(i));
0617 rHeights.push_back(proj->GetBinContent(i));
0618 }
0619 }
0620
0621 if (rPeaks.size() < 5)
0622 {
0623 if (side)
0624 {
0625 m_reco_RPeaks[1].clear();
0626 }
0627 else
0628 {
0629 m_reco_RPeaks[0].clear();
0630 }
0631
0632 if (Verbosity())
0633 {
0634 std::cout << (side ? "side 1" : "side 0") << " has fewer than 5 radial peaks (only has " << rPeaks.size() << "). Returning empty vectors" << std::endl;
0635 }
0636 return hitMatches;
0637 }
0638
0639 double threshold = 0.75;
0640
0641 for (int i = 0; i < (int) rPeaks.size(); i++)
0642 {
0643 std::vector<int> tmpR;
0644 tmpR.push_back(i);
0645
0646 bool closePeak = false;
0647 int currentPeak = -1;
0648
0649 if (rPeaks[i] > 41.0)
0650 {
0651 threshold = 1.0;
0652 }
0653
0654 for (int j = 0; j < (int) finalRPeaks.size(); j++)
0655 {
0656 for (int k = 0; k < (int) groupR[j].size(); k++)
0657 {
0658 if (fabs(rPeaks[i] - rPeaks[groupR[j][k]]) <= threshold || fabs(rPeaks[i] - finalRPeaks[j]) <= threshold)
0659 {
0660 closePeak = true;
0661 currentPeak = j;
0662 break;
0663 }
0664 }
0665 if (closePeak)
0666 {
0667 break;
0668 }
0669 }
0670
0671 if (!closePeak)
0672 {
0673 finalRPeaks.push_back(rPeaks[i]);
0674 finalRHeights.push_back(rHeights[i]);
0675 groupR.push_back(tmpR);
0676 tmpR.clear();
0677 continue;
0678 }
0679
0680 groupR[currentPeak].push_back(i);
0681 double num = 0.0;
0682 double den = 0.0;
0683 for (int j : groupR[currentPeak])
0684 {
0685 double rHeight = proj->GetBinContent(proj->FindBin(rPeaks[j]));
0686 num += rPeaks[j] * rHeight;
0687 den += rHeight;
0688 }
0689
0690 finalRPeaks[currentPeak] = num / den;
0691 finalRHeights[currentPeak] = den;
0692 }
0693
0694 if (side)
0695 {
0696 m_reco_RPeaks[1].clear();
0697
0698 for (double& finalRPeak : finalRPeaks)
0699 {
0700 m_reco_RPeaks[1].push_back(finalRPeak);
0701 }
0702 }
0703 else
0704 {
0705 m_reco_RPeaks[0].clear();
0706
0707 for (double& finalRPeak : finalRPeaks)
0708 {
0709 m_reco_RPeaks[0].push_back(finalRPeak);
0710 }
0711 }
0712
0713 if (Verbosity())
0714 {
0715 std::cout << "finalRPeaks: {";
0716 for (int i = 0; i < (int) finalRPeaks.size() - 1; i++)
0717 {
0718 std::cout << finalRPeaks[i] << ", ";
0719 }
0720 std::cout << finalRPeaks[finalRPeaks.size() - 1] << "}" << std::endl;
0721
0722 if (side)
0723 {
0724 std::cout << "m_reco_RPeaks[1]: {";
0725 for (int i = 0; i < (int) m_reco_RPeaks[1].size() - 1; i++)
0726 {
0727 std::cout << m_reco_RPeaks[1][i] << ", ";
0728 }
0729 std::cout << m_reco_RPeaks[1][m_reco_RPeaks[1].size() - 1] << "}" << std::endl;
0730 }
0731 else
0732 {
0733 std::cout << "m_reco_RPeaks[0]: {";
0734 for (int i = 0; i < (int) m_reco_RPeaks[0].size() - 1; i++)
0735 {
0736 std::cout << m_reco_RPeaks[0][i] << ", ";
0737 }
0738 std::cout << m_reco_RPeaks[0][m_reco_RPeaks[0].size() - 1] << "}" << std::endl;
0739 }
0740
0741 std::cout << "rHeights: {";
0742 for (int i = 0; i < (int) finalRHeights.size() - 1; i++)
0743 {
0744 std::cout << finalRHeights[i] << ", ";
0745 }
0746 std::cout << finalRHeights[finalRHeights.size() - 1] << "}" << std::endl;
0747 }
0748
0749
0750
0751
0752
0753
0754
0755
0756
0757
0758
0759
0760
0761
0762
0763
0764
0765
0766
0767
0768
0769
0770
0771
0772
0773
0774 double first_NN = 100000000.0;
0775 int first_match = -1;
0776 for (int i = 0; i < (int) m_truth_RPeaks.size(); i++)
0777 {
0778 if (fabs(finalRPeaks[0] - m_truth_RPeaks[i]) < first_NN)
0779 {
0780 first_NN = fabs(finalRPeaks[0] - m_truth_RPeaks[i]);
0781 first_match = i;
0782 }
0783 }
0784
0785 double last_NN = 100000000.0;
0786 int last_match = -1;
0787 for (int i = 0; i < (int) m_truth_RPeaks.size(); i++)
0788 {
0789 if (fabs(finalRPeaks[finalRPeaks.size() - 1] - m_truth_RPeaks[i]) < last_NN)
0790 {
0791 last_NN = fabs(finalRPeaks[finalRPeaks.size() - 1] - m_truth_RPeaks[i]);
0792 last_match = i;
0793 }
0794 }
0795
0796 double bestSum = 100000000.0;
0797 int match_i = 0;
0798 int match_j = 0;
0799
0800 int minI = -3;
0801 if (first_match + minI < 0)
0802 {
0803 minI = -first_match;
0804 }
0805
0806 int maxJ = 3;
0807 if (last_match + maxJ >= (int) m_truth_RPeaks.size())
0808 {
0809 maxJ = ((int) m_truth_RPeaks.size()) - 1 - last_match;
0810 }
0811
0812 std::vector<std::vector<double>> matches(3 - minI + 1, std::vector<double>(maxJ + 3 + 1, -999.0));
0813 std::vector<std::vector<std::vector<int>>> NNMatches(3 - minI + 1, std::vector<std::vector<int>>(maxJ + 3 + 1, std::vector<int>(finalRPeaks.size(), -1)));
0814
0815 for (int i = minI; i <= 3; i++)
0816 {
0817 for (int j = -3; j <= maxJ; j++)
0818 {
0819 bool goodPair = true;
0820
0821 if (m_fixShifts)
0822 {
0823 if (m_fieldOn)
0824 {
0825 if (i != -1 || j != -1)
0826 {
0827 goodPair = false;
0828 }
0829 }
0830
0831 if (!m_fieldOn)
0832 {
0833 if (i != 0 || j != 0)
0834 {
0835 goodPair = false;
0836 }
0837 }
0838 }
0839
0840 if (!goodPair)
0841 {
0842 continue;
0843 }
0844
0845 double sum = 0.0;
0846 double m = (m_truth_RPeaks[last_match + j] - m_truth_RPeaks[first_match + i]) / (finalRPeaks[finalRPeaks.size() - 1] - finalRPeaks[0]);
0847 double b = ((m_truth_RPeaks[first_match + i] * finalRPeaks[finalRPeaks.size() - 1]) - (m_truth_RPeaks[last_match + j] * finalRPeaks[0])) / (finalRPeaks[finalRPeaks.size() - 1] - finalRPeaks[0]);
0848
0849 int recoIndex = 0;
0850 for (double finalRPeak : finalRPeaks)
0851 {
0852 int minMatch = 0;
0853 double minResidual = 1000000000.0;
0854 for (int k = 0; k < (int) m_truth_RPeaks.size(); k++)
0855 {
0856 double residual = fabs((m * finalRPeak + b) - m_truth_RPeaks[k]);
0857 if (residual < minResidual)
0858 {
0859 minResidual = residual;
0860 minMatch = k;
0861 if (minResidual < 0.5)
0862 {
0863 break;
0864 }
0865 }
0866 }
0867 sum += fabs((m * finalRPeak + b) - m_truth_RPeaks[minMatch]);
0868
0869 NNMatches[i - minI][j + 3][recoIndex] = minMatch;
0870 recoIndex++;
0871 }
0872
0873
0874 matches[i - minI][j + 3] = sum;
0875 if (m_savehistograms)
0876 {
0877 if (side)
0878 {
0879 m_matchResiduals[1]->Fill(i, j, sum);
0880 }
0881 else
0882 {
0883 m_matchResiduals[0]->Fill(i, j, sum);
0884 }
0885 }
0886 if (sum < bestSum)
0887 {
0888 bestSum = sum;
0889 match_i = i - minI;
0890 match_j = j + 3;
0891 }
0892 }
0893 }
0894
0895 if (Verbosity())
0896 {
0897 std::cout << "best total residual = " << bestSum << " at i " << match_i + minI << " j " << match_j - 3 << std::endl;
0898 for (int i = 0; i < (int) matches.size(); i++)
0899 {
0900 for (int j = 0; j < (int) matches[i].size(); j++)
0901 {
0902 if (matches[i][j] == -999.0)
0903 {
0904 continue;
0905 }
0906 std::cout << "total Residual match i=" << i + minI << " j=" << j - 3 << " : " << matches[i][j] << std::endl;
0907 }
0908 }
0909 }
0910
0911 for (int& i : NNMatches[match_i][match_j])
0912 {
0913 hitMatches.push_back(i);
0914 }
0915
0916 if (Verbosity())
0917 {
0918 for (int i = 0; i < (int) NNMatches.size(); i++)
0919 {
0920 for (int j = 0; j < (int) NNMatches[i].size(); j++)
0921 {
0922 if (matches[i][j] == -999.0)
0923 {
0924 continue;
0925 }
0926 double m = (m_truth_RPeaks[last_match + j - 3] - m_truth_RPeaks[first_match + i + minI]) / (finalRPeaks[finalRPeaks.size() - 1] - finalRPeaks[0]);
0927 double b = ((m_truth_RPeaks[first_match + i + minI] * finalRPeaks[finalRPeaks.size() - 1]) - (m_truth_RPeaks[last_match + j - 3] * finalRPeaks[0])) / (finalRPeaks[finalRPeaks.size() - 1] - finalRPeaks[0]);
0928 for (int k = 0; k < (int) NNMatches[i][j].size(); k++)
0929 {
0930 std::cout << "i " << i + minI << " j " << j - 3 << " Reco index " << k << " recoR=" << finalRPeaks[k] << " shifted R=" << (m * finalRPeaks[k] + b) << " matchIndex=" << NNMatches[i][j][k] << " truth R=" << m_truth_RPeaks[NNMatches[i][j][k]] << " residual=" << (m * finalRPeaks[k] + b) - m_truth_RPeaks[NNMatches[i][j][k]] << std::endl;
0931 }
0932 }
0933 }
0934 }
0935
0936 double final_m = (m_truth_RPeaks[last_match + match_j - 3] - m_truth_RPeaks[first_match + match_i + minI]) / (finalRPeaks[finalRPeaks.size() - 1] - finalRPeaks[0]);
0937 double final_b = ((m_truth_RPeaks[first_match + match_i + minI] * finalRPeaks[finalRPeaks.size() - 1]) - (m_truth_RPeaks[last_match + match_j - 3] * finalRPeaks[0])) / (finalRPeaks[finalRPeaks.size() - 1] - finalRPeaks[0]);
0938
0939 if (side)
0940 {
0941 m_m[1] = final_m;
0942 m_b[1] = final_b;
0943 m_matchLow[1] = match_i + minI;
0944 m_matchHigh[1] = match_j - 3;
0945 }
0946 else
0947 {
0948 m_m[0] = final_m;
0949 m_b[0] = final_b;
0950 m_matchLow[0] = match_i + minI;
0951 m_matchHigh[0] = match_j - 3;
0952 }
0953
0954 return hitMatches;
0955 }
0956
0957 int TpcCentralMembraneMatching::getClusterRMatch(double clusterR, int side)
0958 {
0959 double closestDist = 100.;
0960 int closestPeak = -1;
0961
0962
0963 for (int j = 0; j < (int) m_reco_RPeaks[side].size(); j++)
0964 {
0965 if (std::abs(clusterR - m_reco_RPeaks[side][j]) < closestDist)
0966 {
0967 closestDist = std::abs(clusterR - m_reco_RPeaks[side][j]);
0968 closestPeak = j;
0969 }
0970 }
0971
0972
0973 if (closestPeak != -1)
0974 {
0975 return m_reco_RMatches[side][closestPeak];
0976 }
0977 else
0978 {
0979 return -1;
0980 }
0981 }
0982
0983
0984 int TpcCentralMembraneMatching::InitRun(PHCompositeNode* topNode)
0985 {
0986 if(!m_fieldOn)
0987 {
0988 m_useHeader = false;
0989 }
0990
0991 if (m_savehistograms)
0992 {
0993 static constexpr float max_dr = 5.0;
0994 static constexpr float max_dphi = 0.05;
0995
0996 fout = new TFile(m_histogramfilename.c_str(), "RECREATE");
0997
0998 hxy_reco = new TH2F("hxy_reco", "reco cluster x:y", 800, -100, 100, 800, -80, 80);
0999 hxy_truth = new TH2F("hxy_truth", "truth cluster x:y", 800, -100, 100, 800, -80, 80);
1000
1001 hdrdphi = new TH2F("hdrdphi", "dr vs dphi", 800, -max_dr, max_dr, 800, -max_dphi, max_dphi);
1002 hdrdphi->GetXaxis()->SetTitle("dr");
1003 hdrdphi->GetYaxis()->SetTitle("dphi");
1004
1005 hrdr = new TH2F("hrdr", "dr vs r", 800, 0.0, 80.0, 800, -max_dr, max_dr);
1006 hrdr->GetXaxis()->SetTitle("r");
1007 hrdr->GetYaxis()->SetTitle("dr");
1008
1009 hrdphi = new TH2F("hrdphi", "dphi vs r", 800, 0.0, 80.0, 800, -max_dphi, max_dphi);
1010 hrdphi->GetXaxis()->SetTitle("r");
1011 hrdphi->GetYaxis()->SetTitle("dphi");
1012
1013 hdphi = new TH1F("hdphi", "dph", 800, -max_dphi, max_dphi);
1014 hdphi->GetXaxis()->SetTitle("dphi");
1015
1016 hdr1_single = new TH1F("hdr1_single", "innner dr single", 200, -max_dr, max_dr);
1017 hdr2_single = new TH1F("hdr2_single", "mid dr single", 200, -max_dr, max_dr);
1018 hdr3_single = new TH1F("hdr3_single", "outer dr single", 200, -max_dr, max_dr);
1019 hdr1_double = new TH1F("hdr1_double", "innner dr double", 200, -max_dr, max_dr);
1020 hdr2_double = new TH1F("hdr2_double", "mid dr double", 200, -max_dr, max_dr);
1021 hdr3_double = new TH1F("hdr3_double", "outer dr double", 200, -max_dr, max_dr);
1022 hdrphi = new TH1F("hdrphi", "r * dphi", 200, -0.05, 0.05);
1023 hnclus = new TH1F("hnclus", " nclusters ", 3, 0., 3.);
1024
1025 m_debugfile = new TFile(m_debugfilename.c_str(), "RECREATE");
1026
1027 match_tree = new TTree("match_tree", "Match TTree");
1028
1029 match_tree->Branch("event", &m_event_index);
1030 match_tree->Branch("matched", &m_matched);
1031 match_tree->Branch("truthIndex", &m_truthIndex);
1032 match_tree->Branch("truthR", &m_truthR);
1033 match_tree->Branch("truthPhi", &m_truthPhi);
1034 match_tree->Branch("recoR", &m_recoR);
1035 match_tree->Branch("recoPhi", &m_recoPhi);
1036 match_tree->Branch("recoZ", &m_recoZ);
1037 match_tree->Branch("rawR", &m_rawR);
1038 match_tree->Branch("rawPhi", &m_rawPhi);
1039 match_tree->Branch("staticR", &m_staticR);
1040 match_tree->Branch("staticPhi", &m_staticPhi);
1041 match_tree->Branch("fitMode", &m_fitMode);
1042 match_tree->Branch("side", &m_side);
1043 match_tree->Branch("adc", &m_adc);
1044 match_tree->Branch("nhits", &m_nhits);
1045 match_tree->Branch("nLayers", &m_nLayers);
1046 match_tree->Branch("nIPhi", &m_nIPhi);
1047 match_tree->Branch("nIT", &m_nIT);
1048 match_tree->Branch("layerSD", &m_layersSD);
1049 match_tree->Branch("IPhiSD", &m_IPhiSD);
1050 match_tree->Branch("ITSD", &m_ITSD);
1051 match_tree->Branch("layerWeightedSD", &m_layersWeightedSD);
1052 match_tree->Branch("IPhiWeightedSD", &m_IPhiWeightedSD);
1053 match_tree->Branch("ITWeightedSD", &m_ITWeightedSD);
1054 match_tree->Branch("lowShift", &m_lowShift);
1055 match_tree->Branch("highShift", &m_highShift);
1056 match_tree->Branch("phiRotation", &m_phiRotation);
1057 match_tree->Branch("distanceToTruth", &m_distanceToTruth);
1058 match_tree->Branch("NNDistance", &m_NNDistance);
1059 match_tree->Branch("NNR", &m_NNR);
1060 match_tree->Branch("NNPhi", &m_NNPhi);
1061 match_tree->Branch("NNIndex", &m_NNIndex);
1062
1063
1064 m_matchResiduals[0] = new TH2F("matchResiduals_0", "Matching Residuals TPC South;Shift of smallest R from NN match;Shift of largest R from NN match", 7, -3.5, 3.5, 7, -3.5, 3.5);
1065 m_matchResiduals[1] = new TH2F("matchResiduals_1", "Matching Residuals TPC North;Shift of smallest R from NN match;Shift of largest R from NN match", 7, -3.5, 3.5, 7, -3.5, 3.5);
1066 }
1067
1068 truth_r_phi[0] = new TH2F("truth_r_phi_0", "truth r vs #phi side 0;#phi (rad); r (cm)", 360, -M_PI, M_PI, 500, 0, 100);
1069 truth_r_phi[1] = new TH2F("truth_r_phi_1", "truth r vs #phi side 1;#phi (rad); r (cm)", 360, -M_PI, M_PI, 500, 0, 100);
1070
1071 reco_r_phi[0] = new TH2F("reco_r_phi_0", "reco R vs #phi side 0;#phi (rad); r (cm)", 360, -M_PI, M_PI, 350, 20, 90);
1072 reco_r_phi[1] = new TH2F("reco_r_phi_1", "reco R vs #phi side 1;#phi (rad); r (cm)", 360, -M_PI, M_PI, 350, 20, 90);
1073
1074
1075
1076
1077 const double phi_petal = M_PI / 9.0;
1078
1079
1080
1081
1082
1083
1084
1085 auto save_truth_position = [&](TVector3 source)
1086 {
1087 source.SetZ(-1);
1088 source.RotateZ(-M_PI / 18);
1089 m_truth_pos.push_back(source);
1090 truth_r_phi[0]->Fill(source.Phi(), source.Perp());
1091
1092 source.SetZ(+1);
1093 source.RotateZ(M_PI / 18);
1094 source.SetX(-1*source.X());
1095 m_truth_pos.push_back(source);
1096 truth_r_phi[1]->Fill(source.Phi(), source.Perp());
1097
1098 };
1099
1100
1101 for (int j = 0; j < nRadii; ++j)
1102 {
1103 for (int i = 0; i < nGoodStripes_R1_e[j]; ++i)
1104 {
1105 for (int k = 0; k < 18; ++k)
1106 {
1107 TVector3 dummyPos(cx1_e[i][j], cy1_e[i][j], 0.0);
1108 dummyPos.RotateZ(k * phi_petal);
1109 save_truth_position(dummyPos);
1110
1111
1112 int truth_index_1;
1113 int truth_index_0;
1114 if(k<=8)
1115 {
1116 truth_index_0 = (26-k)*10000 + j*100 + (nGoodStripes_R1_e[j]-i-1);
1117 truth_index_1 = (8-k)*10000 + j*100 + (nGoodStripes_R1_e[j]-i-1);
1118 }
1119 else
1120 {
1121 truth_index_0 = (44-k)*10000 + j*100 + (nGoodStripes_R1_e[j]-i-1);
1122 truth_index_1 = (26-k)*10000 + j*100 + (nGoodStripes_R1_e[j]-i-1);
1123 }
1124 m_truth_index.push_back(truth_index_0);
1125 m_truth_index.push_back(truth_index_1);
1126
1127 if (Verbosity() > 2)
1128 {
1129 std::cout << " i " << i << " j " << j << " k " << k << " x1 " << dummyPos.X() << " y1 " << dummyPos.Y()
1130 << " theta " << std::atan2(dummyPos.Y(), dummyPos.X())
1131 << " radius " << get_r(dummyPos.X(), dummyPos.y()) << std::endl;
1132 }
1133 if (m_savehistograms)
1134 {
1135 hxy_truth->Fill(dummyPos.X(), dummyPos.Y());
1136 }
1137 }
1138 }
1139 }
1140
1141
1142 for (int j = 0; j < nRadii; ++j)
1143 {
1144 for (int i = 0; i < nGoodStripes_R1[j]; ++i)
1145 {
1146 for (int k = 0; k < 18; ++k)
1147 {
1148 TVector3 dummyPos(cx1[i][j], cy1[i][j], 0.0);
1149 dummyPos.RotateZ(k * phi_petal);
1150 save_truth_position(dummyPos);
1151
1152
1153 int truth_index_1;
1154 int truth_index_0;
1155 if(k<=8)
1156 {
1157 truth_index_0 = (26-k)*10000 + (j+8)*100 + (nGoodStripes_R1[j]-i-1);
1158 truth_index_1 = (8-k)*10000 + (j+8)*100 + (nGoodStripes_R1[j]-i-1);
1159 }
1160 else
1161 {
1162 truth_index_0 = (44-k)*10000 + (j+8)*100 + (nGoodStripes_R1[j]-i-1);
1163 truth_index_1 = (26-k)*10000 + (j+8)*100 + (nGoodStripes_R1[j]-i-1);
1164 }
1165 m_truth_index.push_back(truth_index_0);
1166 m_truth_index.push_back(truth_index_1);
1167
1168 if (Verbosity() > 2)
1169 {
1170 std::cout << " i " << i << " j " << j << " k " << k << " x1 " << dummyPos.X() << " y1 " << dummyPos.Y()
1171 << " theta " << std::atan2(dummyPos.Y(), dummyPos.X())
1172 << " radius " << get_r(dummyPos.X(), dummyPos.y()) << std::endl;
1173 }
1174 if (m_savehistograms)
1175 {
1176 hxy_truth->Fill(dummyPos.X(), dummyPos.Y());
1177 }
1178 }
1179 }
1180 }
1181
1182 for (int j = 0; j < nRadii; ++j)
1183 {
1184 for (int i = 0; i < nGoodStripes_R2[j]; ++i)
1185 {
1186 for (int k = 0; k < 18; ++k)
1187 {
1188 TVector3 dummyPos(cx2[i][j], cy2[i][j], 0.0);
1189 dummyPos.RotateZ(k * phi_petal);
1190 save_truth_position(dummyPos);
1191
1192
1193 int truth_index_1;
1194 int truth_index_0;
1195 if(k<=8)
1196 {
1197 truth_index_0 = (26-k)*10000 + (j+16)*100 + (nGoodStripes_R2[j]-i-1);
1198 truth_index_1 = (8-k)*10000 + (j+16)*100 + (nGoodStripes_R2[j]-i-1);
1199 }
1200 else
1201 {
1202 truth_index_0 = (44-k)*10000 + (j+16)*100 + (nGoodStripes_R2[j]-i-1);
1203 truth_index_1 = (26-k)*10000 + (j+16)*100 + (nGoodStripes_R2[j]-i-1);
1204 }
1205 m_truth_index.push_back(truth_index_0);
1206 m_truth_index.push_back(truth_index_1);
1207
1208 if (Verbosity() > 2)
1209 {
1210 std::cout << " i " << i << " j " << j << " k " << k << " x1 " << dummyPos.X() << " y1 " << dummyPos.Y()
1211 << " theta " << std::atan2(dummyPos.Y(), dummyPos.X())
1212 << " radius " << get_r(dummyPos.X(), dummyPos.y()) << std::endl;
1213 }
1214 if (m_savehistograms)
1215 {
1216 hxy_truth->Fill(dummyPos.X(), dummyPos.Y());
1217 }
1218 }
1219 }
1220 }
1221
1222 for (int j = 0; j < nRadii; ++j)
1223 {
1224 for (int i = 0; i < nGoodStripes_R3[j]; ++i)
1225 {
1226 for (int k = 0; k < 18; ++k)
1227 {
1228 TVector3 dummyPos(cx3[i][j], cy3[i][j], 0.0);
1229 dummyPos.RotateZ(k * phi_petal);
1230 save_truth_position(dummyPos);
1231
1232
1233 int truth_index_1;
1234 int truth_index_0;
1235 if(k<=8)
1236 {
1237 truth_index_0 = (26-k)*10000 + (j+24)*100 + (nGoodStripes_R3[j]-i-1);
1238 truth_index_1 = (8-k)*10000 + (j+24)*100 + (nGoodStripes_R3[j]-i-1);
1239 }
1240 else
1241 {
1242 truth_index_0 = (44-k)*10000 + (j+24)*100 + (nGoodStripes_R3[j]-i-1);
1243 truth_index_1 = (26-k)*10000 + (j+24)*100 + (nGoodStripes_R3[j]-i-1);
1244 }
1245 m_truth_index.push_back(truth_index_0);
1246 m_truth_index.push_back(truth_index_1);
1247
1248 if (Verbosity() > 2)
1249 {
1250 std::cout << " i " << i << " j " << j << " k " << k << " x1 " << dummyPos.X() << " y1 " << dummyPos.Y()
1251 << " theta " << std::atan2(dummyPos.Y(), dummyPos.X())
1252 << " radius " << get_r(dummyPos.X(), dummyPos.y()) << std::endl;
1253 }
1254 if (m_savehistograms)
1255 {
1256 hxy_truth->Fill(dummyPos.X(), dummyPos.Y());
1257 }
1258 }
1259 }
1260 }
1261
1262 int ret = GetNodes(topNode);
1263 return ret;
1264 }
1265
1266
1267 int TpcCentralMembraneMatching::process_event(PHCompositeNode* topNode)
1268 {
1269 std::vector<TVector3> reco_pos;
1270 std::vector<TVector3> static_pos;
1271 std::vector<TVector3> raw_pos;
1272 std::vector<bool> reco_side;
1273 std::vector<bool> fitMode;
1274 std::vector<unsigned int> reco_nhits;
1275 std::vector<unsigned int> reco_adc;
1276 std::vector<unsigned int> reco_nLayers;
1277 std::vector<unsigned int> reco_nIPhi;
1278 std::vector<unsigned int> reco_nIT;
1279 std::vector<float> reco_SDLayer;
1280 std::vector<float> reco_SDIPhi;
1281 std::vector<float> reco_SDIT;
1282 std::vector<float> reco_SDWeightedLayer;
1283 std::vector<float> reco_SDWeightedIPhi;
1284 std::vector<float> reco_SDWeightedIT;
1285
1286 eventHeader = findNode::getClass<EventHeader>(topNode, "EventHeader");
1287 if(!eventHeader)
1288 {
1289 std::cout << PHWHERE << " EventHeader Node missing, abort" << std::endl;
1290 return Fun4AllReturnCodes::ABORTRUN;
1291 }
1292
1293 if(m_useHeader && eventHeader->get_EvtSequence() == 0)
1294 {
1295 m_useHeader = false;
1296 }
1297
1298
1299 if(m_useHeader)
1300 {
1301 m_event_index = eventHeader->get_EvtSequence();
1302 }
1303
1304 if (!m_corrected_CMcluster_map || m_corrected_CMcluster_map->size() < 1000)
1305 {
1306 if(!m_useHeader)
1307 {
1308 m_event_index++;
1309 }
1310 return Fun4AllReturnCodes::EVENT_OK;
1311 }
1312
1313
1314 if (Verbosity())
1315 {
1316 std::cout << PHWHERE << " working on event " << m_event_index << " which has " << m_corrected_CMcluster_map->size() << " clusters" << std::endl;
1317 }
1318
1319
1320 for (const auto& harray : {m_dcc_out->m_hDRint, m_dcc_out->m_hDPint, m_dcc_out->m_hDZint, m_dcc_out->m_hentries})
1321 {
1322 for (const auto& h : harray)
1323 {
1324 h->Reset();
1325 }
1326 }
1327
1328 reco_r_phi[0]->Reset();
1329 reco_r_phi[1]->Reset();
1330
1331 int nClus_gtMin = 0;
1332 int clusterIndex = 0;
1333
1334
1335 auto clusrange = m_corrected_CMcluster_map->getClusters();
1336 for (auto cmitr = clusrange.first;
1337 cmitr != clusrange.second;
1338 ++cmitr)
1339 {
1340 const auto& [cmkey, cmclus_orig] = *cmitr;
1341
1342 LaserCluster* cmclus = cmclus_orig;
1343 const unsigned int nhits = cmclus->getNhits();
1344
1345 const unsigned int adc = cmclus->getAdc();
1346 bool side = (bool) TpcDefs::getSide(cmkey);
1347
1348
1349
1350
1351
1352
1353 nClus_gtMin++;
1354
1355
1356
1357
1358 Acts::Vector3 pos(cmclus->getX(), cmclus->getY(), (side ? 1.0 : -1.0));
1359 TVector3 tmp_raw(pos[0], pos[1], pos[2]);
1360 if (m_dcc_in_module_edge)
1361 {
1362 pos = m_distortionCorrection.get_corrected_position(pos, m_dcc_in_module_edge);
1363 }
1364 if (m_dcc_in_static)
1365 {
1366 pos = m_distortionCorrection.get_corrected_position(pos, m_dcc_in_static);
1367 }
1368 TVector3 tmp_static(pos[0], pos[1], pos[2]);
1369
1370 if (m_dcc_in_average)
1371 {
1372 pos = m_distortionCorrection.get_corrected_position(pos, m_dcc_in_average);
1373 }
1374 TVector3 tmp_pos(pos[0], pos[1], pos[2]);
1375
1376
1377
1378
1379 reco_pos.push_back(tmp_pos);
1380 static_pos.push_back(tmp_static);
1381 raw_pos.push_back(tmp_raw);
1382 reco_side.push_back(side);
1383 fitMode.push_back(cmclus->getFitMode());
1384 reco_nhits.push_back(nhits);
1385 reco_adc.push_back(adc);
1386 reco_nLayers.push_back(cmclus->getNLayers());
1387 reco_nIPhi.push_back(cmclus->getNIPhi());
1388 reco_nIT.push_back(cmclus->getNIT());
1389 reco_SDLayer.push_back(cmclus->getSDLayer());
1390 reco_SDIPhi.push_back(cmclus->getSDIPhi());
1391 reco_SDIT.push_back(cmclus->getSDIT());
1392 reco_SDWeightedLayer.push_back(cmclus->getSDWeightedLayer());
1393 reco_SDWeightedIPhi.push_back(cmclus->getSDWeightedIPhi());
1394 reco_SDWeightedIT.push_back(cmclus->getSDWeightedIT());
1395
1396 if (side == 0)
1397 {
1398 reco_r_phi[0]->Fill(tmp_pos.Phi(), tmp_pos.Perp());
1399 }
1400 else
1401 {
1402 reco_r_phi[1]->Fill(tmp_pos.Phi(), tmp_pos.Perp());
1403 }
1404
1405 if (Verbosity())
1406 {
1407 double raw_rad = sqrt(cmclus->getX() * cmclus->getX() + cmclus->getY() * cmclus->getY());
1408 double static_rad = sqrt(tmp_static.X() * tmp_static.X() + tmp_static.Y() * tmp_static.Y());
1409 double corr_rad = sqrt(tmp_pos.X() * tmp_pos.X() + tmp_pos.Y() * tmp_pos.Y());
1410 std::cout << "cluster " << clusterIndex << std::endl;
1411 clusterIndex++;
1412 std::cout << "found raw cluster " << cmkey << " side " << side << " with x " << cmclus->getX() << " y " << cmclus->getY() << " z " << cmclus->getZ() << " radius " << raw_rad << std::endl;
1413 std::cout << " --- static corrected positions: " << tmp_static.X() << " " << tmp_static.Y() << " " << tmp_static.Z() << " radius " << static_rad << std::endl;
1414 std::cout << " --- corrected positions: " << tmp_pos.X() << " " << tmp_pos.Y() << " " << tmp_pos.Z() << " radius " << corr_rad << std::endl;
1415 }
1416
1417 if (m_savehistograms)
1418 {
1419 hxy_reco->Fill(tmp_pos.X(), tmp_pos.Y());
1420 }
1421 }
1422
1423 if (nClus_gtMin < 25)
1424 {
1425 if(!m_useHeader)
1426 {
1427 m_event_index++;
1428 }
1429 return Fun4AllReturnCodes::EVENT_OK;
1430 }
1431
1432 int truth_index = 0;
1433 int nMatched = 0;
1434 std::vector<bool> truth_matched(m_truth_pos.size(), false);
1435 std::vector<bool> reco_matched(reco_pos.size(), false);
1436 std::vector<int> truth_matchedRecoIndex(m_truth_pos.size(), -1);
1437 std::vector<int> reco_matchedTruthIndex(reco_pos.size(), -1);
1438
1439 std::vector<float> reco_distToNN(reco_pos.size(), 10000.0);
1440 std::vector<float> NNDist(reco_pos.size(), 100000.0);
1441 std::vector<float> NNR(reco_pos.size(), 100000.0);
1442 std::vector<float> NNPhi(reco_pos.size(), 100000.0);
1443 std::vector<int> NNIndex(reco_pos.size(), -1);
1444
1445 std::vector<std::vector<int>> truth_NNRecoIndex(m_truth_pos.size(), std::vector<int>(0));
1446
1447
1448 if(m_doFancy){
1449
1450
1451 m_recoRotation[0][1] = getPhiRotation_smoothed(truth_r_phi[0]->ProjectionX("hR2", 206, 290), reco_r_phi[0]->ProjectionX("cR2_0", 206, 290), false);
1452 m_recoRotation[0][0] = getPhiRotation_smoothed(truth_r_phi[0]->ProjectionX("hR1", 151, 206), reco_r_phi[0]->ProjectionX("cR1_0", 151, 206), false);
1453 m_recoRotation[0][2] = getPhiRotation_smoothed(truth_r_phi[0]->ProjectionX("hR3", 290, 499), reco_r_phi[0]->ProjectionX("cR3_0", 290, 499), false);
1454
1455 m_recoRotation[1][1] = getPhiRotation_smoothed(truth_r_phi[1]->ProjectionX("hR2", 206, 290), reco_r_phi[1]->ProjectionX("cR2_1", 206, 290), true);
1456 m_recoRotation[1][0] = getPhiRotation_smoothed(truth_r_phi[1]->ProjectionX("hR1", 151, 206), reco_r_phi[1]->ProjectionX("cR1_1", 151, 206), true);
1457 m_recoRotation[1][2] = getPhiRotation_smoothed(truth_r_phi[1]->ProjectionX("hR3", 290, 499), reco_r_phi[1]->ProjectionX("cR3_1", 290, 499), true);
1458
1459 m_reco_RMatches[0] = doGlobalRMatching(reco_r_phi[0], false);
1460 m_reco_RMatches[1] = doGlobalRMatching(reco_r_phi[1], true);
1461
1462 if (Verbosity())
1463 {
1464 for (int i = 0; i < (int) m_reco_RMatches[0].size(); i++)
1465 {
1466 std::cout << "side 0 cluster index " << i << " hit match " << m_reco_RMatches[0][i] << " recoPeak=" << m_reco_RPeaks[0][i] << " shifted recoPeak=" << (m_m[0] * m_reco_RPeaks[0][i] + m_b[0]) << " truthPeak=" << m_truth_RPeaks[m_reco_RMatches[0][i]] << " residual=" << m_truth_RPeaks[m_reco_RMatches[0][i]] - (m_m[0] * m_reco_RPeaks[0][i] + m_b[0]) << std::endl;
1467 }
1468
1469 for (int i = 0; i < (int) m_reco_RMatches[1].size(); i++)
1470 {
1471 std::cout << "side 1 cluster index " << i << " hit match " << m_reco_RMatches[1][i] << " recoPeak=" << m_reco_RPeaks[1][i] << " shifted recoPeak=" << (m_m[1] * m_reco_RPeaks[1][i] + m_b[1]) << " truthPeak=" << m_truth_RPeaks[m_reco_RMatches[1][i]] << " residual=" << m_truth_RPeaks[m_reco_RMatches[1][i]] - (m_m[1] * m_reco_RPeaks[1][i] + m_b[1]) << std::endl;
1472 }
1473 }
1474
1475
1476
1477
1478 for (const auto& truth : m_truth_pos)
1479 {
1480 double tR = get_r(truth.X(), truth.Y());
1481 double tPhi = truth.Phi();
1482 double tZ = truth.Z();
1483
1484 int truthRIndex = -1;
1485
1486
1487 for (int k = 0; k < (int) m_truth_RPeaks.size(); k++)
1488 {
1489 if (std::abs(tR - m_truth_RPeaks[k]) < 0.5)
1490 {
1491 truthRIndex = k;
1492 break;
1493 }
1494 }
1495
1496 if (truthRIndex == -1)
1497 {
1498 truth_index++;
1499 continue;
1500 }
1501
1502 double prev_dphi = 10000.0;
1503
1504 int recoMatchIndex = -1;
1505 unsigned int reco_index = 0;
1506 for (const auto& reco : reco_pos)
1507 {
1508 if (reco_matched[reco_index] || reco_nhits[reco_index] < m_nHitsInCuster_minimum)
1509 {
1510 reco_index++;
1511 continue;
1512 }
1513
1514 double rR = get_r(reco.X(), reco.Y());
1515 double rPhi = reco.Phi();
1516 bool side = reco_side[reco_index];
1517
1518 int region = -1;
1519
1520 if (rR < 41)
1521 {
1522 region = 0;
1523 }
1524 else if (rR >= 41 && rR < 58)
1525 {
1526 region = 1;
1527 }
1528 else if (rR >= 58)
1529 {
1530 region = 2;
1531 }
1532
1533 if (region != -1)
1534 {
1535 if (side)
1536 {
1537 rPhi -= m_recoRotation[1][region];
1538 }
1539 else
1540 {
1541 rPhi -= m_recoRotation[0][region];
1542 }
1543 }
1544
1545 int clustRMatchIndex = -1;
1546 clustRMatchIndex = getClusterRMatch(rR, (int) (side ? 1 : 0));
1547
1548 if (clustRMatchIndex == -1 || truthRIndex != clustRMatchIndex)
1549 {
1550 reco_index++;
1551 continue;
1552 }
1553
1554 if ((!side && tZ > 0) || (side && tZ < 0))
1555 {
1556 reco_index++;
1557 continue;
1558 }
1559
1560 auto dphi = delta_phi(tPhi - rPhi);
1561 if (fabs(dphi) > m_phi_cut)
1562 {
1563 reco_index++;
1564 continue;
1565 }
1566
1567 if (fabs(dphi) < fabs(prev_dphi))
1568 {
1569 prev_dphi = dphi;
1570 recoMatchIndex = reco_index;
1571 truth_matched[truth_index] = true;
1572 }
1573
1574 reco_index++;
1575 }
1576
1577 if (recoMatchIndex != -1)
1578 {
1579 truth_matchedRecoIndex[truth_index] = recoMatchIndex;
1580 reco_matched[recoMatchIndex] = true;
1581 reco_matchedTruthIndex[recoMatchIndex] = truth_index;
1582 nMatched++;
1583
1584 if (Verbosity() > 2)
1585 {
1586 std::cout << "truth " << truth_index << " matched to reco " << recoMatchIndex << " and there are now " << nMatched << " matches" << std::endl;
1587 std::cout << "tR=" << std::setw(10) << tR << " tPhi=" << std::setw(10) << tPhi << " tZ=" << std::setw(10) << tZ << std::endl;
1588 std::cout << "rR=" << std::setw(10) << get_r(reco_pos[recoMatchIndex].X(), reco_pos[recoMatchIndex].Y()) << " rPhi=" << std::setw(10) << reco_pos[recoMatchIndex].Phi() << " rZ=" << std::setw(10) << reco_pos[recoMatchIndex].Z() << " rSide=" << reco_side[recoMatchIndex] << " dPhi=" << prev_dphi << std::endl;
1589 }
1590 }
1591
1592 truth_index++;
1593 }
1594
1595
1596
1597 int recoIndex = 0;
1598 for (const auto& reco : reco_pos)
1599 {
1600 double rR = get_r(reco.X(), reco.Y());
1601 double rPhi = reco.Phi();
1602 bool side = reco_side[recoIndex];
1603
1604 int region = -1;
1605
1606 if (rR < 41)
1607 {
1608 region = 0;
1609 }
1610 else if (rR >= 41 && rR < 58)
1611 {
1612 region = 1;
1613 }
1614 else if (rR >= 58)
1615 {
1616 region = 2;
1617 }
1618
1619 if (region != -1)
1620 {
1621 if (side)
1622 {
1623 rPhi -= m_recoRotation[1][region];
1624 }
1625 else
1626 {
1627 rPhi -= m_recoRotation[0][region];
1628 }
1629 }
1630
1631 int clustRMatchIndex = -1;
1632 clustRMatchIndex = getClusterRMatch(rR, (int) (side ? 1 : 0));
1633
1634
1635 int truthMatchIndex = -1;
1636 truth_index = 0;
1637 for (const auto& truth : m_truth_pos)
1638 {
1639 double tR = get_r(truth.X(), truth.Y());
1640 double tPhi = truth.Phi();
1641 double tZ = truth.Z();
1642
1643 if( (!side && tZ > 0) || (side && tZ < 0) )
1644 {
1645 truth_index++;
1646 continue;
1647 }
1648
1649 if(sqrt(pow(truth.X() - reco.X(),2) + pow(truth.Y() - reco.Y(),2)) < NNDist[recoIndex])
1650 {
1651 NNDist[recoIndex] = sqrt(pow(truth.X() - reco.X(),2) + pow(truth.Y() - reco.Y(),2));
1652 NNR[recoIndex] = tR;
1653 NNPhi[recoIndex] = tPhi;
1654 NNIndex[recoIndex] = m_truth_index[truth_index];
1655 }
1656
1657 if (clustRMatchIndex == -1)
1658 {
1659 continue;
1660 }
1661
1662 if (reco_matched[recoIndex])
1663 {
1664 if (reco_matchedTruthIndex[recoIndex] == truth_index)
1665 {
1666 reco_distToNN[recoIndex] = sqrt(pow(truth.X() - reco.X(), 2) + pow(truth.Y() - reco.Y(), 2));
1667 break;
1668 }
1669 else
1670 {
1671 truth_index++;
1672 continue;
1673 }
1674 }
1675
1676 int truthRIndex = -1;
1677
1678
1679 for (int k = 0; k < (int) m_truth_RPeaks.size(); k++)
1680 {
1681 if (std::abs(tR - m_truth_RPeaks[k]) < 0.5)
1682 {
1683 truthRIndex = k;
1684 break;
1685 }
1686 }
1687
1688 if (truthRIndex == -1 || truthRIndex != clustRMatchIndex)
1689 {
1690 truth_index++;
1691 continue;
1692 }
1693
1694 if ((!side && tZ > 0) || (side && tZ < 0))
1695 {
1696 truth_index++;
1697 continue;
1698 }
1699
1700 auto dphi = delta_phi(tPhi - rPhi);
1701 if (fabs(dphi) > m_phi_cut)
1702 {
1703 truth_index++;
1704 continue;
1705 }
1706
1707 float dist = sqrt(pow(tR * sin(tPhi) - rR * sin(rPhi), 2) + pow(tR * cos(tPhi) - rR * cos(rPhi), 2));
1708
1709 if (dist < reco_distToNN[recoIndex])
1710 {
1711 reco_distToNN[recoIndex] = dist;
1712 truthMatchIndex = truth_index;
1713 }
1714
1715 truth_index++;
1716
1717 }
1718
1719 if (!reco_matched[recoIndex])
1720 {
1721 reco_matchedTruthIndex[recoIndex] = truthMatchIndex;
1722 }
1723 recoIndex++;
1724 }
1725 }
1726 else
1727 {
1728 int reco_index = 0;
1729 for (const auto& reco : reco_pos)
1730 {
1731 double rR = get_r(reco.X(), reco.Y());
1732 double rPhi = reco.Phi();
1733 bool side = reco_side[reco_index];
1734
1735 truth_index = 0;
1736 double minNNDist = 100000.0;
1737 int match_localTruth = -1;
1738 for (const auto& truth : m_truth_pos)
1739 {
1740
1741 double tR = get_r(truth.X(), truth.Y());
1742 double tPhi = truth.Phi();
1743 double tZ = truth.Z();
1744
1745 if ((!side && tZ > 0) || (side && tZ < 0))
1746 {
1747 truth_index++;
1748 continue;
1749 }
1750
1751 auto dR = fabs(tR - rR);
1752 if (dR > 5.0)
1753 {
1754 truth_index++;
1755 continue;
1756 }
1757
1758 auto dphi = delta_phi(tPhi - rPhi);
1759 if (fabs(dphi) > 0.05)
1760 {
1761 truth_index++;
1762 continue;
1763 }
1764
1765 double dist = sqrt(pow(truth.X() - reco.X(),2) + pow(truth.Y() - reco.Y(),2));
1766 if (dist < minNNDist)
1767 {
1768 minNNDist = dist;
1769 match_localTruth = truth_index;
1770 }
1771 truth_index++;
1772 }
1773
1774 if(match_localTruth == -1)
1775 {
1776 reco_index++;
1777 continue;
1778 }
1779
1780 truth_NNRecoIndex[match_localTruth].push_back(reco_index);
1781 NNDist[reco_index] = sqrt(pow(m_truth_pos[match_localTruth].X() - reco.X(),2) + pow(m_truth_pos[match_localTruth].Y() - reco.Y(),2));
1782 NNR[reco_index] = get_r(m_truth_pos[match_localTruth].X(), m_truth_pos[match_localTruth].Y());
1783 NNPhi[reco_index] = m_truth_pos[match_localTruth].Phi();
1784 NNIndex[reco_index] = m_truth_index[match_localTruth];
1785
1786 reco_index++;
1787 }
1788
1789 truth_index = 0;
1790 for(const auto& truthIndex : truth_NNRecoIndex)
1791 {
1792 double maxADC = 0.0;
1793 int recoMatchIndex = -1;
1794 for(const auto& recoIndex : truthIndex)
1795 {
1796
1797 if(reco_nLayers[recoIndex] < 2 || reco_nLayers[recoIndex] > 3) continue;
1798
1799 int lowRMatch = 0;
1800 int highRMatch = (int)m_truth_RPeaks.size()-1;
1801 int upperBound = (int)m_truth_RPeaks.size();
1802 while(lowRMatch <= highRMatch)
1803 {
1804 int mid = lowRMatch + (highRMatch - lowRMatch)/2;
1805 if(m_truth_RPeaks[mid] >= reco_pos[recoIndex].Perp())
1806 {
1807 upperBound = mid;
1808 highRMatch = mid - 1;
1809 }
1810 else
1811 {
1812 lowRMatch = mid + 1;
1813 }
1814 }
1815 if(upperBound < 1) upperBound = 1;
1816
1817 double dR = m_truth_pos[truth_index].Perp() - reco_pos[recoIndex].Perp();
1818 if(fabs(dR) > 0.5 * (m_truth_RPeaks[upperBound] - m_truth_RPeaks[upperBound-1]))
1819 {
1820 continue;
1821 }
1822
1823 double dphi = delta_phi(m_truth_pos[truth_index].Phi() - reco_pos[recoIndex].Phi());
1824 if(fabs(dphi) > 0.5*phiSpacing[upperBound-1])
1825 {
1826 continue;
1827 }
1828
1829 if(reco_adc[recoIndex] > maxADC)
1830 {
1831 maxADC = reco_adc[recoIndex];
1832 recoMatchIndex = recoIndex;
1833 }
1834 }
1835
1836 if(recoMatchIndex >= 0)
1837 {
1838 truth_matchedRecoIndex[truth_index] = recoMatchIndex;
1839 reco_matched[recoMatchIndex] = true;
1840 reco_matchedTruthIndex[recoMatchIndex] = truth_index;
1841 truth_matched[truth_index] = true;
1842 nMatched++;
1843 }
1844
1845 truth_index++;
1846 }
1847 }
1848
1849
1850 if (Verbosity())
1851 {
1852 const auto n_valid_truth = std::count_if(m_truth_pos.begin(), m_truth_pos.end(), [](const TVector3& pos)
1853 { return get_r(pos.x(), pos.y()) > 30; });
1854 std::cout << "TpcCentralMembraneMatching::process_event - m_truth_pos size: " << m_truth_pos.size() << std::endl;
1855 std::cout << "TpcCentralMembraneMatching::process_event - m_truth_pos size, r>30cm: " << n_valid_truth << std::endl;
1856 std::cout << "TpcCentralMembraneMatching::process_event - reco_pos size: " << reco_pos.size() << std::endl;
1857 std::cout << "TpcCentralMembraneMatching::process_event - matched_pair size: " << nMatched << std::endl;
1858 }
1859
1860 if (m_savehistograms)
1861 {
1862 for (int i = 0; i < (int) reco_pos.size(); i++)
1863 {
1864 m_matched = reco_matched[i];
1865 m_truthR = m_truth_pos[reco_matchedTruthIndex[i]].Perp();
1866 m_truthPhi = m_truth_pos[reco_matchedTruthIndex[i]].Phi();
1867 m_distanceToTruth = sqrt(pow(m_truth_pos[reco_matchedTruthIndex[i]].Y() - reco_pos[i].Y(), 2) + pow(m_truth_pos[reco_matchedTruthIndex[i]].X() - reco_pos[i].X(), 2));
1868 m_truthIndex = m_truth_index[reco_matchedTruthIndex[i]];
1869 m_recoR = reco_pos[i].Perp();
1870 m_recoPhi = reco_pos[i].Phi();
1871 m_recoZ = reco_pos[i].Z();
1872 m_rawR = raw_pos[i].Perp();
1873 m_rawPhi = raw_pos[i].Phi();
1874 m_staticR = static_pos[i].Perp();
1875 m_staticPhi = static_pos[i].Phi();
1876 m_side = reco_side[i];
1877 m_fitMode = fitMode[i];
1878 m_adc = reco_adc[i];
1879 m_nhits = reco_nhits[i];
1880 m_nLayers = reco_nLayers[i];
1881 m_nIPhi = reco_nIPhi[i];
1882 m_nIT = reco_nIT[i];
1883 m_layersSD = reco_SDLayer[i];
1884 m_IPhiSD = reco_SDIPhi[i];
1885 m_ITSD = reco_SDIT[i];
1886 m_layersWeightedSD = reco_SDWeightedLayer[i];
1887 m_IPhiWeightedSD = reco_SDWeightedIPhi[i];
1888 m_ITWeightedSD = reco_SDWeightedIT[i];
1889
1890
1891 if (m_side)
1892 {
1893 m_lowShift = m_matchLow[1];
1894 m_highShift = m_matchHigh[1];
1895 }
1896 else
1897 {
1898 m_lowShift = m_matchLow[0];
1899 m_highShift = m_matchHigh[0];
1900 }
1901
1902 int region = -1;
1903
1904 if (m_recoR < 41)
1905 {
1906 region = 0;
1907 }
1908 else if (m_recoR >= 41 && m_recoR < 58)
1909 {
1910 region = 1;
1911 }
1912 else if (m_recoR >= 58)
1913 {
1914 region = 2;
1915 }
1916
1917 if (region == -1)
1918 {
1919 m_phiRotation = -999.0;
1920 }
1921 else
1922 {
1923 m_phiRotation = m_recoRotation[m_side][region];
1924 }
1925
1926 m_NNDistance = NNDist[i];
1927 m_NNR = NNR[i];
1928 m_NNPhi = NNPhi[i];
1929 m_NNIndex = NNIndex[i];
1930
1931 match_tree->Fill();
1932 }
1933 }
1934
1935 unsigned int ckey = 0;
1936 for (unsigned int i = 0; i < m_truth_pos.size(); i++)
1937 {
1938 if (!truth_matched[i])
1939 {
1940 continue;
1941 }
1942
1943
1944
1945 int reco_index = truth_matchedRecoIndex[i];
1946
1947 if (reco_index == -1)
1948 {
1949 continue;
1950 }
1951
1952
1953
1954 auto cmdiff = new CMFlashDifferencev1();
1955 cmdiff->setTruthPhi(m_truth_pos[i].Phi());
1956 cmdiff->setTruthR(m_truth_pos[i].Perp());
1957 cmdiff->setTruthZ(m_truth_pos[i].Z());
1958
1959 if (m_averageMode)
1960 {
1961 cmdiff->setRecoPhi(static_pos[reco_index].Phi());
1962 cmdiff->setRecoR(static_pos[reco_index].Perp());
1963 cmdiff->setRecoZ(static_pos[reco_index].Z());
1964 cmdiff->setNclusters(reco_nhits[reco_index]);
1965 }
1966 else
1967 {
1968 cmdiff->setRecoPhi(reco_pos[reco_index].Phi());
1969 cmdiff->setRecoR(reco_pos[reco_index].Perp());
1970 cmdiff->setRecoZ(reco_pos[reco_index].Z());
1971 cmdiff->setNclusters(reco_nhits[reco_index]);
1972 }
1973
1974 if (Verbosity() > 1)
1975 {
1976 std::cout << "adding cmdiff to container with key " << ckey << " in event " << m_event_index << std::endl;
1977 }
1978
1979 m_cm_flash_diffs->addDifferenceSpecifyKey(ckey, cmdiff);
1980
1981
1982 const double clus_r = reco_pos[reco_index].Perp();
1983 double clus_phi = reco_pos[reco_index].Phi();
1984 if (clus_phi < 0)
1985 {
1986 clus_phi += 2 * M_PI;
1987 }
1988
1989
1990 const bool side = reco_side[reco_index];
1991
1992
1993
1994
1995 double dr = reco_pos[reco_index].Perp() - m_truth_pos[i].Perp();
1996 double dphi = delta_phi(reco_pos[reco_index].Phi() - m_truth_pos[i].Phi());
1997 double rdphi = reco_pos[reco_index].Perp() * dphi;
1998 if (m_averageMode)
1999 {
2000 dr = static_pos[reco_index].Perp() - m_truth_pos[i].Perp();
2001 dphi = delta_phi(static_pos[reco_index].Phi() - m_truth_pos[i].Phi());
2002 rdphi = static_pos[reco_index].Perp() * dphi;
2003 }
2004
2005
2006
2007 const double dz = 0.0;
2008
2009
2010
2011
2012
2013
2014
2015 for (const auto& dcc : {m_dcc_out, m_dcc_out_aggregated.get()})
2016 {
2017 static_cast<TH2*>(dcc->m_hDRint[side])->Fill(clus_phi, clus_r, dr);
2018 static_cast<TH2*>(dcc->m_hDPint[side])->Fill(clus_phi, clus_r, rdphi);
2019 static_cast<TH2*>(dcc->m_hDZint[side])->Fill(clus_phi, clus_r, dz);
2020 static_cast<TH2*>(dcc->m_hentries[side])->Fill(clus_phi, clus_r);
2021 }
2022
2023 ckey++;
2024 }
2025
2026 if (Verbosity())
2027 {
2028 std::cout << "TpcCentralMembraneMatching::process_events - cmclusters: " << m_corrected_CMcluster_map->size() << std::endl;
2029 std::cout << "TpcCentralMembraneMatching::process_events - matched pairs: " << nMatched << std::endl;
2030 std::cout << "TpcCentralMembraneMatching::process_events - differences: " << m_cm_flash_diffs->size() << std::endl;
2031 std::cout << "TpcCentralMembraneMatching::process_events - entries: " << m_dcc_out->m_hentries[0]->GetEntries() << ", " << m_dcc_out->m_hentries[1]->GetEntries() << std::endl;
2032 }
2033
2034
2035 normalize_distortions(m_dcc_out);
2036 fill_guarding_bins(m_dcc_out);
2037
2038 if (Verbosity() > 2)
2039 {
2040
2041 auto diffrange = m_cm_flash_diffs->getDifferences();
2042 for (auto cmitr = diffrange.first;
2043 cmitr != diffrange.second;
2044 ++cmitr)
2045 {
2046 auto key = cmitr->first;
2047 auto cmreco = cmitr->second;
2048
2049 std::cout << " key " << key
2050 << " nclus " << cmreco->getNclusters()
2051 << " truth Phi " << cmreco->getTruthPhi() << " reco Phi " << cmreco->getRecoPhi()
2052 << " truth R " << cmreco->getTruthR() << " reco R " << cmreco->getRecoR()
2053 << " truth Z " << cmreco->getTruthZ() << " reco Z " << cmreco->getRecoZ()
2054 << std::endl;
2055 }
2056 }
2057
2058 if(!m_useHeader){
2059 m_event_index++;
2060 }
2061
2062 return Fun4AllReturnCodes::EVENT_OK;
2063 }
2064
2065
2066 int TpcCentralMembraneMatching::End(PHCompositeNode* )
2067 {
2068
2069 if (m_dcc_out_aggregated)
2070 {
2071
2072 if(m_doHadd)
2073 {
2074 normalize_distortions(m_dcc_out_aggregated.get());
2075 }
2076 fill_guarding_bins(m_dcc_out_aggregated.get());
2077
2078
2079 std::unique_ptr<TFile> outputfile(TFile::Open(m_outputfile.c_str(), "RECREATE"));
2080 outputfile->cd();
2081
2082
2083 for (unsigned int i = 0; i < 2; ++i)
2084 {
2085 for (const auto& h : {m_dcc_out_aggregated->m_hDRint[i], m_dcc_out_aggregated->m_hDPint[i], m_dcc_out_aggregated->m_hDZint[i], m_dcc_out_aggregated->m_hentries[i]})
2086 {
2087 if (h)
2088 {
2089 std::cout << "writing in " << m_outputfile << " histogram " << h->GetName() << std::endl;
2090 h->Write();
2091 }
2092 }
2093 }
2094 }
2095
2096
2097 if (m_savehistograms && fout)
2098 {
2099 fout->cd();
2100
2101 hxy_reco->Write();
2102 hxy_truth->Write();
2103 hdrdphi->Write();
2104 hrdr->Write();
2105 hrdphi->Write();
2106 hdphi->Write();
2107 hdrphi->Write();
2108 hdr1_single->Write();
2109 hdr2_single->Write();
2110 hdr3_single->Write();
2111 hdr1_double->Write();
2112 hdr2_double->Write();
2113 hdr3_double->Write();
2114 hnclus->Write();
2115
2116 fout->Close();
2117 }
2118
2119 if (m_savehistograms && m_debugfile)
2120 {
2121 m_debugfile->cd();
2122
2123 std::cout << "writing in " << m_debugfilename << " match tree " << std::endl;
2124
2125
2126
2127
2128
2129
2130
2131
2132
2133
2134
2135
2136
2137
2138 m_debugfile->Write();
2139 std::cout << "closing " << m_debugfilename << std::endl;
2140 m_debugfile->Close();
2141 std::cout << m_debugfilename << " has been closed" << std::endl;
2142
2143 }
2144
2145 return Fun4AllReturnCodes::EVENT_OK;
2146 }
2147
2148
2149
2150 int TpcCentralMembraneMatching::GetNodes(PHCompositeNode* topNode)
2151 {
2152
2153
2154
2155
2156
2157 m_corrected_CMcluster_map = findNode::getClass<LaserClusterContainer>(topNode, "LASER_CLUSTER");
2158 if (!m_corrected_CMcluster_map)
2159 {
2160 std::cout << PHWHERE << "CORRECTED_CM_CLUSTER Node missing, abort." << std::endl;
2161 return Fun4AllReturnCodes::ABORTRUN;
2162 }
2163
2164
2165 m_dcc_in_module_edge = findNode::getClass<TpcDistortionCorrectionContainer>(topNode, "TpcDistortionCorrectionContainerModuleEdge");
2166 if (m_dcc_in_module_edge)
2167 {
2168 std::cout << "TpcCentralMembraneMatching::GetNodes - found TPC distortion correction container module edge" << std::endl;
2169 }
2170
2171
2172 m_dcc_in_static = findNode::getClass<TpcDistortionCorrectionContainer>(topNode, "TpcDistortionCorrectionContainerStatic");
2173 if (m_dcc_in_static)
2174 {
2175 std::cout << "TpcCentralMembraneMatching::GetNodes - found TPC distortion correction container static" << std::endl;
2176 }
2177
2178
2179 m_dcc_in_average = findNode::getClass<TpcDistortionCorrectionContainer>(topNode, "TpcDistortionCorrectionContainerAverage");
2180 if (m_dcc_in_average)
2181 {
2182 std::cout << "TpcCentralMembraneMatching::GetNodes - found TPC distortion correction container average" << std::endl;
2183 }
2184
2185 PHNodeIterator iter(topNode);
2186
2187
2188 PHCompositeNode* dstNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "DST"));
2189 if (!dstNode)
2190 {
2191 std::cout << PHWHERE << "DST Node missing, doing nothing." << std::endl;
2192 return Fun4AllReturnCodes::ABORTRUN;
2193 }
2194
2195
2196
2197
2198
2199
2200
2201
2202
2203
2204
2205
2206 auto flashDiffContainer = findNode::getClass<CMFlashDifferenceContainerv1>(topNode, "CM_FLASH_DIFFERENCES");
2207 if (!flashDiffContainer)
2208 {
2209 PHNodeIterator dstiter(dstNode);
2210 PHCompositeNode* DetNode =
2211 dynamic_cast<PHCompositeNode*>(dstiter.findFirst("PHCompositeNode", "TRKR"));
2212 if (!DetNode)
2213 {
2214 DetNode = new PHCompositeNode("TRKR");
2215 dstNode->addNode(DetNode);
2216 }
2217
2218 flashDiffContainer = new CMFlashDifferenceContainerv1;
2219 PHIODataNode<PHObject>* CMFlashDifferenceNode =
2220 new PHIODataNode<PHObject>(flashDiffContainer, "CM_FLASH_DIFFERENCES", "PHObject");
2221
2222 DetNode->addNode(CMFlashDifferenceNode);
2223 }
2224
2225
2226 m_cm_flash_diffs = findNode::getClass<CMFlashDifferenceContainerv1>(topNode, "CM_FLASH_DIFFERENCES");
2227 if (!m_cm_flash_diffs)
2228 {
2229 std::cout << PHWHERE << " ERROR: Can't find CM_FLASH_DIFFERENCES." << std::endl;
2230 return Fun4AllReturnCodes::ABORTRUN;
2231 }
2232
2233
2234
2235
2236
2237
2238
2239
2240
2241
2242
2243
2244
2245
2246
2247
2248
2249
2250
2251
2252
2253
2254
2255
2256
2257
2258 const std::string dcc_out_node_name = "TpcDistortionCorrectionContainerFluctuation";
2259 m_dcc_out = findNode::getClass<TpcDistortionCorrectionContainer>(topNode, dcc_out_node_name);
2260 if (!m_dcc_out)
2261 {
2262
2263 auto runNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "RUN"));
2264 if (!runNode)
2265 {
2266 std::cout << "TpcCentralMembraneMatching::InitRun - RUN Node missing, quitting" << std::endl;
2267 return Fun4AllReturnCodes::ABORTRUN;
2268 }
2269
2270 std::cout << "TpcCentralMembraneMatching::GetNodes - creating TpcDistortionCorrectionContainer in node " << dcc_out_node_name << std::endl;
2271 m_dcc_out = new TpcDistortionCorrectionContainer;
2272 m_dcc_out->m_dimensions = 2;
2273 m_dcc_out->m_phi_hist_in_radians = false;
2274 m_dcc_out->m_interpolate_z = true;
2275 auto node = new PHDataNode<TpcDistortionCorrectionContainer>(m_dcc_out, dcc_out_node_name);
2276 runNode->addNode(node);
2277 }
2278
2279
2280
2281
2282
2283 m_dcc_out_aggregated.reset(new TpcDistortionCorrectionContainer);
2284
2285
2286 const float phiMin = m_phiMin - (m_phiMax - m_phiMin) / m_phibins;
2287 const float phiMax = m_phiMax + (m_phiMax - m_phiMin) / m_phibins;
2288
2289 const float rMin = m_rMin - (m_rMax - m_rMin) / m_rbins;
2290 const float rMax = m_rMax + (m_rMax - m_rMin) / m_rbins;
2291
2292
2293
2294
2295
2296
2297
2298
2299
2300
2301
2302
2303
2304
2305
2306
2307
2308
2309
2310
2311
2312
2313
2314
2315
2316
2317
2318
2319
2320
2321
2322
2323
2324
2325 const std::array<const std::string, 2> extension = {{"_negz", "_posz"}};
2326 for (const auto& dcc : {m_dcc_out, m_dcc_out_aggregated.get()})
2327 {
2328
2329 dcc->m_dimensions = 2;
2330
2331
2332 for (int i = 0; i < 2; ++i)
2333 {
2334 delete dcc->m_hDPint[i];
2335 dcc->m_hDPint[i] = new TH2F((boost::format("hIntDistortionP%s") % extension[i]).str().c_str(), (boost::format("hIntDistortionP%s") % extension[i]).str().c_str(), m_phibins + 2, phiMin, phiMax, m_rbins + 2, rMin, rMax);
2336 delete dcc->m_hDRint[i];
2337 dcc->m_hDRint[i] = new TH2F((boost::format("hIntDistortionR%s") % extension[i]).str().c_str(), (boost::format("hIntDistortionR%s") % extension[i]).str().c_str(), m_phibins + 2, phiMin, phiMax, m_rbins + 2, rMin, rMax);
2338 delete dcc->m_hDZint[i];
2339 dcc->m_hDZint[i] = new TH2F((boost::format("hIntDistortionZ%s") % extension[i]).str().c_str(), (boost::format("hIntDistortionZ%s") % extension[i]).str().c_str(), m_phibins + 2, phiMin, phiMax, m_rbins + 2, rMin, rMax);
2340 delete dcc->m_hentries[i];
2341 dcc->m_hentries[i] = new TH2I((boost::format("hEntries%s") % extension[i]).str().c_str(), (boost::format("hEntries%s") % extension[i]).str().c_str(), m_phibins + 2, phiMin, phiMax, m_rbins + 2, rMin, rMax);
2342
2343
2344
2345
2346
2347 }
2348 }
2349
2350 return Fun4AllReturnCodes::EVENT_OK;
2351 }
2352
2353
2354 void TpcCentralMembraneMatching::CalculateCenters(
2355 int nPads,
2356 const std::array<double, nRadii>& R,
2357 std::array<int, nRadii>& nGoodStripes,
2358 const std::array<int, nRadii>& keepUntil,
2359 std::array<int, nRadii>& nStripesIn,
2360 std::array<int, nRadii>& nStripesBefore,
2361 double cx[][nRadii], double cy[][nRadii])
2362 {
2363 const double phi_module = M_PI / 6.0;
2364 const int pr_mult = 3;
2365 const int dw_mult = 8;
2366 const double diffwidth = 0.6 * mm;
2367 const double adjust = 0.015;
2368
2369 double theta = 0.0;
2370
2371
2372
2373
2374 std::array<double, nRadii> spacing{};
2375 for (int i = 0; i < nRadii; i++)
2376 {
2377 spacing[i] = 2.0 * ((dw_mult * diffwidth / R[i]) + (pr_mult * phi_module / nPads));
2378 }
2379
2380
2381 for (int j = 0; j < nRadii; j++)
2382 {
2383 int i_out = 0;
2384 for (int i = keepThisAndAfter[j]; i < keepUntil[j]; i++)
2385 {
2386 if (j % 2 == 0)
2387 {
2388 theta = i * spacing[j] + (spacing[j] / 2) - adjust;
2389 cx[i_out][j] = R[j] * cos(theta) / cm;
2390 cy[i_out][j] = R[j] * sin(theta) / cm;
2391 }
2392 else
2393 {
2394 theta = (i + 1) * spacing[j] - adjust;
2395 cx[i_out][j] = R[j] * cos(theta) / cm;
2396 cy[i_out][j] = R[j] * sin(theta) / cm;
2397 }
2398
2399 if (Verbosity() > 2)
2400 {
2401 std::cout << " j " << j << " i " << i << " i_out " << i_out << " theta " << theta << " cx " << cx[i_out][j] << " cy " << cy[i_out][j]
2402 << " radius " << sqrt(pow(cx[i_out][j], 2) + pow(cy[i_out][j], 2)) << std::endl;
2403 }
2404
2405 i_out++;
2406
2407 nStripesBefore_R1_e[0] = 0;
2408
2409 nStripesIn[j] = keepUntil[j] - keepThisAndAfter[j];
2410 if (j == 0)
2411 {
2412 nStripesBefore[j] = 0;
2413 }
2414 else
2415 {
2416 nStripesBefore[j] = nStripesIn[j - 1] + nStripesBefore[j - 1];
2417 }
2418 nStripesBefore_R1_e[0] = 0;
2419 }
2420 nGoodStripes[j] = i_out;
2421 }
2422 }