File indexing completed on 2025-08-06 08:19:30
0001 #include "PHG4TpcDirectLaser.h"
0002
0003 #include <phparameter/PHParameterInterface.h> // for PHParameterInterface
0004
0005 #include <g4main/PHG4HitContainer.h>
0006 #include <g4main/PHG4HitDefs.h> // for get_volume_id
0007 #include <g4main/PHG4Hitv1.h>
0008 #include <g4main/PHG4Particlev3.h>
0009 #include <g4main/PHG4TruthInfoContainer.h>
0010 #include <g4main/PHG4VtxPointv1.h>
0011
0012 #include <fun4all/Fun4AllReturnCodes.h>
0013 #include <fun4all/SubsysReco.h> // for SubsysReco
0014
0015 #include <phool/PHCompositeNode.h>
0016 #include <phool/PHIODataNode.h>
0017 #include <phool/PHNode.h>
0018 #include <phool/PHNodeIterator.h>
0019 #include <phool/PHObject.h> // for PHObject
0020 #include <phool/getClass.h>
0021 #include <phool/phool.h> // for PHWHERE
0022
0023 #include <trackbase_historic/SvtxTrackMap.h>
0024 #include <trackbase_historic/SvtxTrackMap_v2.h>
0025 #include <trackbase_historic/SvtxTrack_v2.h>
0026
0027 #include <TFile.h>
0028 #include <TNtuple.h>
0029 #include <TVector3.h> // for TVector3, operator*
0030
0031 #include <gsl/gsl_const_mksa.h> // for the speed of light
0032
0033 #include <cassert>
0034 #include <iostream> // for operator<<, basic_os...
0035 #include <optional>
0036
0037 namespace
0038 {
0039 using PHG4Particle_t = PHG4Particlev3;
0040 using PHG4VtxPoint_t = PHG4VtxPointv1;
0041 using PHG4Hit_t = PHG4Hitv1;
0042
0043
0044 template <class T>
0045 inline constexpr T square(const T& x)
0046 {
0047 return x * x;
0048 }
0049
0050
0051 static const int detId = PHG4HitDefs::get_volume_id("PHG4TpcDirectLaser");
0052
0053
0054
0055 static constexpr double cm = 1.0;
0056
0057
0058
0059 static constexpr double speed_of_light = GSL_CONST_MKSA_SPEED_OF_LIGHT * 1e-7;
0060
0061
0062 static constexpr double maxHitLength = 1. * cm;
0063
0064
0065 static constexpr double halflength_tpc = 105.5 * cm;
0066
0067
0068 static constexpr double begin_CM = 20. * cm;
0069 static constexpr double end_CM = 78. * cm;
0070
0071
0072 static constexpr double halfwidth_CM = 0.5 * cm;
0073
0074
0075 std::optional<TVector3> central_membrane_intersection(const TVector3& start, const TVector3& direction)
0076 {
0077 const double end = start.z() > 0 ? halfwidth_CM : -halfwidth_CM;
0078 const double dist = end - start.z();
0079
0080
0081 if (direction.z() == 0)
0082 {
0083 return std::nullopt;
0084 }
0085
0086
0087 if (dist * direction.z() < 0)
0088 {
0089 return std::nullopt;
0090 }
0091
0092 const double direction_scale = dist / direction.z();
0093 return start + direction * direction_scale;
0094 }
0095
0096
0097 std::optional<TVector3> endcap_intersection(const TVector3& start, const TVector3& direction)
0098 {
0099 const double end = start.z() > 0 ? halflength_tpc : -halflength_tpc;
0100 const double dist = end - start.z();
0101
0102
0103 if (direction.z() == 0)
0104 {
0105 return std::nullopt;
0106 }
0107
0108
0109 if (dist * direction.z() < 0)
0110 {
0111 return std::nullopt;
0112 }
0113
0114 const double direction_scale = dist / direction.z();
0115 return start + direction * direction_scale;
0116 }
0117
0118
0119 std::optional<TVector3> cylinder_line_intersection(const TVector3& s, const TVector3& v, double radius)
0120 {
0121 const double R2 = square(radius);
0122
0123
0124
0125 const double a = square(v.x()) + square(v.y());
0126 const double b = 2 * (v.x() * s.x() + v.y() * s.y());
0127 const double c = square(s.x()) + square(s.y()) - R2;
0128
0129 const double rootterm = square(b) - 4 * a * c;
0130
0131
0132
0133
0134
0135
0136 if (rootterm < 0 || a == 0)
0137 {
0138 return std::nullopt;
0139 }
0140
0141
0142 const double sqrtterm = std::sqrt(rootterm);
0143 const double t1 = (-b + sqrtterm) / (2 * a);
0144 const double t2 = (-b - sqrtterm) / (2 * a);
0145
0146
0147
0148
0149
0150 const double& min_t = (t2 < t1 && t2 > 0) ? t2 : t1;
0151 return s + v * min_t;
0152 }
0153
0154
0155 std::optional<TVector3> field_cage_intersection(const TVector3& start, const TVector3& direction)
0156 {
0157 auto ofc_strike = cylinder_line_intersection(start, direction, end_CM);
0158 auto ifc_strike = cylinder_line_intersection(start, direction, begin_CM);
0159
0160
0161 if (!ifc_strike)
0162 {
0163 return ofc_strike;
0164 }
0165 if (!ofc_strike)
0166 {
0167 return ifc_strike;
0168 }
0169
0170
0171 const auto ifc_dist = (ifc_strike->Z() - start.Z()) / direction.Z();
0172 const auto ofc_dist = (ofc_strike->Z() - start.Z()) / direction.Z();
0173
0174 if (ifc_dist < 0)
0175 {
0176 return (ofc_dist > 0) ? ofc_strike : std::nullopt;
0177 }
0178 else if (ofc_dist < 0)
0179 {
0180 return ifc_strike;
0181 }
0182 else
0183 {
0184 return (ifc_dist < ofc_dist) ? ifc_strike : ofc_strike;
0185 }
0186 }
0187
0188
0189 inline std::ostream& operator<<(std::ostream& out, const TVector3& vector)
0190 {
0191 out << "( " << vector.x() << ", " << vector.y() << ", " << vector.z() << ")";
0192 return out;
0193 }
0194
0195 }
0196
0197
0198 PHG4TpcDirectLaser::PHG4TpcDirectLaser(const std::string& name)
0199 : SubsysReco(name)
0200 , PHParameterInterface(name)
0201 {
0202 InitializeParameters();
0203 }
0204
0205
0206 int PHG4TpcDirectLaser::InitRun(PHCompositeNode* topNode)
0207 {
0208
0209 m_g4truthinfo = findNode::getClass<PHG4TruthInfoContainer>(topNode, "G4TruthInfo");
0210 if (!m_g4truthinfo)
0211 {
0212 std::cout << "Fun4AllDstPileupMerger::load_nodes - creating node G4TruthInfo" << std::endl;
0213
0214 PHNodeIterator iter(topNode);
0215 auto dstNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "DST"));
0216 if (!dstNode)
0217 {
0218 std::cout << PHWHERE << "DST Node missing, aborting." << std::endl;
0219 return Fun4AllReturnCodes::ABORTRUN;
0220 }
0221
0222 m_g4truthinfo = new PHG4TruthInfoContainer();
0223 dstNode->addNode(new PHIODataNode<PHObject>(m_g4truthinfo, "G4TruthInfo", "PHObject"));
0224 }
0225
0226
0227 hitnodename = "G4HIT_" + detector;
0228 auto* g4hit = findNode::getClass<PHG4HitContainer>(topNode, hitnodename.c_str());
0229 if (!g4hit)
0230 {
0231 std::cout << Name() << " Could not locate G4HIT node " << hitnodename << std::endl;
0232 return Fun4AllReturnCodes::ABORTRUN;
0233 }
0234
0235
0236
0237 m_track_map = findNode::getClass<SvtxTrackMap>(topNode, m_track_map_name);
0238 if (!m_track_map)
0239 {
0240
0241 PHNodeIterator iter(topNode);
0242 auto dstNode = static_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "DST"));
0243 if (!dstNode)
0244 {
0245 std::cout << PHWHERE << "DST Node missing, aborting." << std::endl;
0246 return Fun4AllReturnCodes::ABORTRUN;
0247 }
0248
0249
0250 iter = PHNodeIterator(dstNode);
0251 auto node = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "SVTX"));
0252 if (!node)
0253 {
0254 dstNode->addNode(node = new PHCompositeNode("SVTX"));
0255 }
0256
0257
0258 m_track_map = new SvtxTrackMap_v2;
0259 node->addNode(new PHIODataNode<PHObject>(m_track_map, m_track_map_name, "PHObject"));
0260 }
0261
0262
0263 UpdateParametersWithMacro();
0264 electrons_per_cm = get_int_param("electrons_per_cm");
0265 electrons_per_gev = get_double_param("electrons_per_gev");
0266
0267
0268 SetupLasers();
0269
0270
0271 if (m_steppingpattern == true)
0272 {
0273 std::cout << "PHG4TpcDirectLaser::InitRun - m_steppingpattern: " << m_steppingpattern << std::endl;
0274 std::cout << "PHG4TpcDirectLaser::InitRun - nTotalSteps: " << nTotalSteps << std::endl;
0275 }
0276 else
0277 {
0278 std::cout << "PHG4TpcDirectLaser::InitRun - m_autoAdvanceDirectLaser: " << m_autoAdvanceDirectLaser << std::endl;
0279 std::cout << "PHG4TpcDirectLaser::InitRun - phi steps: " << nPhiSteps << " min: " << minPhi << " max: " << maxPhi << std::endl;
0280 std::cout << "PHG4TpcDirectLaser::InitRun - theta steps: " << nThetaSteps << " min: " << minTheta << " max: " << maxTheta << std::endl;
0281 std::cout << "PHG4TpcDirectLaser::InitRun - nTotalSteps: " << nTotalSteps << std::endl;
0282 }
0283 std::cout << "PHG4TpcDirectLaser::InitRun - electrons_per_cm: " << electrons_per_cm << std::endl;
0284 std::cout << "PHG4TpcDirectLaser::InitRun - electrons_per_gev " << electrons_per_gev << std::endl;
0285
0286
0287
0288 std::string LASER_ANGLES_ROOTFILE = std::string(getenv("CALIBRATIONROOT")) + "/TPC/DirectLaser/theta_phi_laser.root";
0289 TFile* infile1 = TFile::Open(LASER_ANGLES_ROOTFILE.c_str());
0290
0291 pattern = (TNtuple*) infile1->Get("angles");
0292 pattern->SetBranchAddress("#theta", &theta_p);
0293 pattern->SetBranchAddress("#phi", &phi_p);
0294
0295 return Fun4AllReturnCodes::EVENT_OK;
0296 }
0297
0298
0299 int PHG4TpcDirectLaser::process_event(PHCompositeNode* topNode)
0300 {
0301
0302 m_g4truthinfo = findNode::getClass<PHG4TruthInfoContainer>(topNode, "G4TruthInfo");
0303 assert(m_g4truthinfo);
0304
0305
0306 m_g4hitcontainer = findNode::getClass<PHG4HitContainer>(topNode, hitnodename.c_str());
0307 assert(m_g4hitcontainer);
0308
0309
0310 m_track_map = findNode::getClass<SvtxTrackMap>(topNode, m_track_map_name);
0311 assert(m_track_map);
0312
0313 if (m_autoAdvanceDirectLaser || m_steppingpattern)
0314 {
0315 AimToNextPatternStep();
0316 }
0317
0318 else
0319 {
0320
0321 AimToThetaPhi(arbitrary_theta, arbitrary_phi);
0322 }
0323
0324 return Fun4AllReturnCodes::EVENT_OK;
0325 }
0326
0327
0328 void PHG4TpcDirectLaser::SetDefaultParameters()
0329 {
0330
0331
0332
0333
0334
0335
0336 static constexpr double Ne_dEdx = 1.56;
0337 static constexpr double CF4_dEdx = 7.00;
0338 static constexpr double Ne_NTotal = 43;
0339 static constexpr double CF4_NTotal = 100;
0340 static constexpr double Tpc_NTot = 0.5 * Ne_NTotal + 0.5 * CF4_NTotal;
0341 static constexpr double Tpc_dEdx = 0.5 * Ne_dEdx + 0.5 * CF4_dEdx;
0342 static constexpr double Tpc_ElectronsPerKeV = Tpc_NTot / Tpc_dEdx;
0343
0344
0345 set_default_double_param("electrons_per_gev", Tpc_ElectronsPerKeV * 1e6);
0346
0347
0348 set_default_int_param("electrons_per_cm", 72);
0349 }
0350
0351
0352 void PHG4TpcDirectLaser::SetPhiStepping(int n, double min, double max)
0353 {
0354 if (n < 0 || max < min)
0355 {
0356 std::cout << PHWHERE << " - invalid" << std::endl;
0357 return;
0358 }
0359 nPhiSteps = n;
0360 minPhi = min;
0361 maxPhi = max;
0362 nTotalSteps = nThetaSteps * nPhiSteps;
0363 return;
0364 }
0365
0366 void PHG4TpcDirectLaser::SetThetaStepping(int n, double min, double max)
0367 {
0368 if (n < 0 || max < min)
0369 {
0370 std::cout << PHWHERE << " - invalid" << std::endl;
0371 return;
0372 }
0373 nThetaSteps = n;
0374 minTheta = min;
0375 maxTheta = max;
0376 nTotalSteps = nThetaSteps * nPhiSteps;
0377
0378 return;
0379 }
0380
0381
0382 void PHG4TpcDirectLaser::SetFileStepping(int n)
0383 {
0384 if (n < 0 || n > 13802)
0385 {
0386 std::cout << PHWHERE << " - invalid" << std::endl;
0387 return;
0388 }
0389 nTotalSteps = n;
0390
0391 return;
0392 }
0393
0394
0395
0396 void PHG4TpcDirectLaser::SetupLasers()
0397 {
0398
0399 m_lasers.clear();
0400
0401
0402 const TVector3 position_base(60 * cm, 0., halflength_tpc);
0403
0404
0405 for (int i = 0; i < 8; ++i)
0406 {
0407 Laser laser;
0408
0409
0410
0411
0412
0413
0414 laser.m_position = position_base;
0415 if (i < 4)
0416 {
0417 laser.m_position.SetZ(position_base.z());
0418 laser.m_direction = -1;
0419 laser.m_phi = M_PI / 2 * i - (15 * M_PI / 180);
0420 }
0421 else
0422 {
0423 laser.m_position.SetZ(-position_base.z());
0424 laser.m_direction = 1;
0425 laser.m_phi = M_PI / 2 * i + (15 * M_PI / 180);
0426 }
0427
0428
0429 laser.m_position.RotateZ(laser.m_phi);
0430
0431
0432 m_lasers.push_back(laser);
0433
0434
0435
0436 }
0437 }
0438
0439
0440 void PHG4TpcDirectLaser::AimToNextPatternStep()
0441 {
0442 if (nTotalSteps >= 1)
0443 {
0444 if (m_steppingpattern)
0445 {
0446 AimToPatternStep_File(currentPatternStep);
0447 ++currentPatternStep;
0448 }
0449 else
0450 {
0451 AimToPatternStep(currentPatternStep);
0452 ++currentPatternStep;
0453 }
0454 }
0455 }
0456
0457
0458 void PHG4TpcDirectLaser::AimToThetaPhi(double theta, double phi)
0459 {
0460 if (Verbosity())
0461 {
0462 std::cout << "PHG4TpcDirectLaser::AimToThetaPhi - theta: " << theta << " phi: " << phi << std::endl;
0463 }
0464
0465
0466 for (const auto& laser : m_lasers)
0467 {
0468 AppendLaserTrack(theta, phi, laser);
0469 }
0470 }
0471
0472
0473 void PHG4TpcDirectLaser::AimToPatternStep(int n)
0474 {
0475
0476 n = n % nTotalSteps;
0477
0478 if (Verbosity())
0479 {
0480 std::cout << "PHG4TpcDirectLaser::AimToPatternStep - step: " << n << "/" << nTotalSteps << std::endl;
0481 }
0482
0483
0484 currentPatternStep = n;
0485
0486
0487 const int thetaStep = n / nPhiSteps;
0488 const double theta = minTheta + thetaStep * (maxTheta - minTheta) / nThetaSteps;
0489
0490
0491 const int phiStep = n % nPhiSteps;
0492 const double phi = minPhi + phiStep * (maxPhi - minPhi) / nPhiSteps;
0493
0494
0495 AimToThetaPhi(theta, phi);
0496
0497 return;
0498 }
0499
0500
0501
0502 void PHG4TpcDirectLaser::AimToPatternStep_File(int n)
0503 {
0504
0505 n = n % nTotalSteps;
0506
0507 if (Verbosity())
0508 {
0509 std::cout << "PHG4TpcDirectLaser::AimToPatternStep_File - step: " << n << "/" << nTotalSteps << std::endl;
0510 }
0511
0512
0513 currentPatternStep = n;
0514
0515 pattern->GetEntry(n);
0516
0517
0518 std::cout << "From file, current entry = " << n << " Theta: " << theta_p << " Phi: " << phi_p << std::endl;
0519
0520 const double theta = theta_p * M_PI / 180.;
0521
0522
0523 const double phi = phi_p * M_PI / 180.;
0524
0525
0526 AimToThetaPhi(theta, phi);
0527
0528 return;
0529 }
0530
0531
0532
0533 void PHG4TpcDirectLaser::AppendLaserTrack(double theta, double phi, const PHG4TpcDirectLaser::Laser& laser)
0534 {
0535 if (!m_g4hitcontainer)
0536 {
0537 std::cout << PHWHERE << "invalid g4hit container. aborting" << std::endl;
0538 return;
0539 }
0540
0541
0542 const auto& pos = laser.m_position;
0543
0544
0545 const auto& direction = laser.m_direction;
0546 TVector3 dir(0, 0, direction);
0547
0548
0549 dir.RotateY(theta * direction);
0550
0551 if (laser.m_direction == -1)
0552 {
0553 dir.RotateZ(phi);
0554 }
0555 else
0556 {
0557 dir.RotateZ(-phi);
0558 }
0559
0560
0561 dir.RotateZ(laser.m_phi);
0562
0563
0564 if (Verbosity())
0565 {
0566 std::cout << "PHG4TpcDirectLaser::AppendLaserTrack - position: " << pos << " direction: " << dir << std::endl;
0567 }
0568
0569
0570 static constexpr double total_momentum = 1;
0571
0572
0573 int trackid = -1;
0574
0575
0576 if (m_g4truthinfo)
0577 {
0578
0579 const auto vtxid = m_g4truthinfo->maxvtxindex() + 1;
0580 const auto vertex = new PHG4VtxPoint_t(pos.x(), pos.y(), pos.z(), 0, vtxid);
0581 m_g4truthinfo->AddVertex(vtxid, vertex);
0582
0583
0584 trackid = m_g4truthinfo->maxtrkindex() + 1;
0585
0586
0587 auto particle = new PHG4Particle_t();
0588 particle->set_track_id(trackid);
0589 particle->set_vtx_id(vtxid);
0590 particle->set_parent_id(0);
0591 particle->set_primary_id(trackid);
0592 particle->set_px(total_momentum * dir.x());
0593 particle->set_py(total_momentum * dir.y());
0594 particle->set_pz(total_momentum * dir.z());
0595
0596 m_g4truthinfo->AddParticle(trackid, particle);
0597 }
0598
0599
0600 if (m_track_map)
0601 {
0602 SvtxTrack_v2 track;
0603 track.set_x(pos.x());
0604 track.set_y(pos.y());
0605 track.set_z(pos.z());
0606
0607
0608 track.set_px(total_momentum * dir.x());
0609 track.set_py(total_momentum * dir.y());
0610 track.set_pz(total_momentum * dir.z());
0611
0612
0613 m_track_map->insert(&track);
0614
0615 if (Verbosity())
0616 {
0617 std::cout << "PHG4TpcDirectLaser::AppendLaserTrack - position: " << pos << " direction: " << dir << std::endl;
0618 }
0619 }
0620
0621
0622
0623
0624
0625
0626
0627 const auto plane_strike = (pos.z() * dir.z() > 0) ? endcap_intersection(pos, dir) : central_membrane_intersection(pos, dir);
0628
0629
0630 const auto fc_strike = field_cage_intersection(pos, dir);
0631
0632
0633 if (!(plane_strike || fc_strike))
0634 {
0635 return;
0636 }
0637
0638
0639
0640 const TVector3& strike = (fc_strike && (!plane_strike || fc_strike->z() / dir.z() < plane_strike->z() / dir.z())) ? *fc_strike : *plane_strike;
0641
0642
0643 TVector3 delta = (strike - pos);
0644 double fullLength = delta.Mag();
0645 int nHitSteps = fullLength / maxHitLength + 1;
0646
0647 TVector3 start = pos;
0648 TVector3 end = start;
0649 TVector3 step = dir * (maxHitLength / (dir.Mag()));
0650 double stepLength = 0;
0651
0652 if (Verbosity())
0653 {
0654 std::cout << "PHG4TpcDirectLaser::AppendLaserTrack -"
0655 << " fullLength: " << fullLength
0656 << " nHitSteps: " << nHitSteps
0657 << std::endl;
0658 }
0659
0660 for (int i = 0; i < nHitSteps; i++)
0661 {
0662 start = end;
0663 if (i + 1 == nHitSteps)
0664 {
0665
0666 end = strike;
0667 delta = start - end;
0668 stepLength = delta.Mag();
0669 }
0670 else
0671 {
0672
0673 end = start + step;
0674 stepLength = step.Mag();
0675 }
0676
0677
0678 auto hit = new PHG4Hit_t;
0679 hit->set_trkid(trackid);
0680 hit->set_layer(99);
0681
0682
0683 hit->set_x(0, start.X() / cm);
0684 hit->set_y(0, start.Y() / cm);
0685 hit->set_z(0, start.Z() / cm);
0686 hit->set_t(0, (start - pos).Mag() / speed_of_light);
0687
0688 hit->set_x(1, end.X() / cm);
0689 hit->set_y(1, end.Y() / cm);
0690 hit->set_z(1, end.Z() / cm);
0691 hit->set_t(1, (end - pos).Mag() / speed_of_light);
0692
0693
0694 hit->set_px(0, dir.X());
0695 hit->set_py(0, dir.Y());
0696 hit->set_pz(0, dir.Z());
0697
0698 hit->set_px(1, dir.X());
0699 hit->set_py(1, dir.Y());
0700 hit->set_pz(1, dir.Z());
0701
0702 const double totalE = electrons_per_cm * stepLength / electrons_per_gev;
0703
0704 hit->set_eion(totalE);
0705 hit->set_edep(totalE);
0706 m_g4hitcontainer->AddHit(detId, hit);
0707 }
0708
0709 return;
0710 }