File indexing completed on 2025-08-05 08:18:12
0001
0002
0003
0004
0005
0006 #include "PHG4MicromegasHitReco.h"
0007
0008 #include <fun4all/Fun4AllReturnCodes.h>
0009 #include <fun4all/SubsysReco.h> // for SubsysReco
0010
0011 #include <g4detectors/PHG4CylinderGeom.h>
0012 #include <g4detectors/PHG4CylinderGeomContainer.h>
0013
0014 #include <g4main/PHG4Hit.h>
0015 #include <g4main/PHG4HitContainer.h>
0016 #include <g4main/PHG4Hitv1.h>
0017
0018 #include <micromegas/CylinderGeomMicromegas.h>
0019 #include <micromegas/MicromegasDefs.h>
0020
0021 #include <phool/PHCompositeNode.h>
0022 #include <phool/PHIODataNode.h>
0023 #include <phool/PHNode.h>
0024 #include <phool/PHNodeIterator.h>
0025 #include <phool/PHObject.h> // for PHObject
0026 #include <phool/PHRandomSeed.h>
0027 #include <phool/getClass.h>
0028 #include <phool/phool.h>
0029
0030 #include <phparameter/PHParameterInterface.h> // for PHParameterInterface
0031
0032 #include <trackbase/ActsGeometry.h>
0033 #include <trackbase/TrkrDefs.h>
0034 #include <trackbase/TrkrHitSet.h>
0035 #include <trackbase/TrkrHitSetContainerv1.h>
0036 #include <trackbase/TrkrHitTruthAssocv1.h>
0037 #include <trackbase/TrkrHitv2.h>
0038
0039 #include <TVector2.h>
0040 #include <TVector3.h>
0041
0042 #include <gsl/gsl_randist.h>
0043 #include <gsl/gsl_rng.h> // for gsl_rng_alloc
0044
0045 #include <cassert>
0046 #include <cmath> // for atan2, sqrt, M_PI
0047 #include <cstdlib> // for exit
0048 #include <iostream> // for operator<<, basic...
0049 #include <map> // for _Rb_tree_const_it...
0050 #include <numeric>
0051
0052 namespace
0053 {
0054
0055 template <class T>
0056 inline constexpr T square(const T& x)
0057 {
0058 return x * x;
0059 }
0060
0061
0062 template <class T>
0063 inline T gaus(const T& x, const T& sigma)
0064 {
0065 return std::exp(-square(x / sigma) / 2) / (sigma * std::sqrt(2 * M_PI));
0066 }
0067
0068
0069 template <class T>
0070 inline T bind_angle(const T& angle)
0071 {
0072 if (angle >= M_PI)
0073 {
0074 return angle - 2 * M_PI;
0075 }
0076 else if (angle < -M_PI)
0077 {
0078 return angle + 2 * M_PI;
0079 }
0080 else
0081 {
0082 return angle;
0083 }
0084 }
0085
0086
0087 template <class T>
0088 inline T get_rectangular_fraction(const T& xloc, const T& sigma, const T& pitch)
0089 {
0090 return (std::erf((xloc + pitch / 2) / (M_SQRT2 * sigma)) - std::erf((xloc - pitch / 2) / (M_SQRT2 * sigma))) / 2;
0091 }
0092
0093
0094
0095
0096
0097 template <class T>
0098 inline T get_zigzag_fraction(const T& xloc, const T& sigma, const T& pitch)
0099 {
0100 return
0101
0102 (pitch - xloc) * (std::erf(xloc / (M_SQRT2 * sigma)) - std::erf((xloc - pitch) / (M_SQRT2 * sigma))) / (pitch * 2) + (gaus(xloc - pitch, sigma) - gaus(xloc, sigma)) * square(sigma) / pitch
0103
0104
0105 + (pitch + xloc) * (std::erf((xloc + pitch) / (M_SQRT2 * sigma)) - std::erf(xloc / (M_SQRT2 * sigma))) / (pitch * 2) + (gaus(xloc + pitch, sigma) - gaus(xloc, sigma)) * square(sigma) / pitch;
0106 }
0107
0108
0109 [[maybe_unused]] inline std::ostream& operator<<(std::ostream& out, const TVector3& position)
0110 {
0111 out << "(" << position.x() << ", " << position.y() << ", " << position.z() << ")";
0112 return out;
0113 }
0114
0115 }
0116
0117
0118 PHG4MicromegasHitReco::PHG4MicromegasHitReco(const std::string& name)
0119 : SubsysReco(name)
0120 , PHParameterInterface(name)
0121 {
0122
0123 const uint seed = PHRandomSeed();
0124 m_rng.reset(gsl_rng_alloc(gsl_rng_mt19937));
0125 gsl_rng_set(m_rng.get(), seed);
0126
0127 InitializeParameters();
0128 }
0129
0130
0131 int PHG4MicromegasHitReco::InitRun(PHCompositeNode* topNode)
0132 {
0133 UpdateParametersWithMacro();
0134
0135
0136 m_tmin = get_double_param("micromegas_tmin");
0137 m_tmax = get_double_param("micromegas_tmax");
0138 m_electrons_per_gev = get_double_param("micromegas_electrons_per_gev");
0139 m_gain = get_double_param("micromegas_gain");
0140 m_cloud_sigma = get_double_param("micromegas_cloud_sigma");
0141 m_diffusion_trans = get_double_param("micromegas_diffusion_trans");
0142 m_added_smear_sigma_z = get_double_param("micromegas_added_smear_sigma_z");
0143 m_added_smear_sigma_rphi = get_double_param("micromegas_added_smear_sigma_rphi");
0144
0145
0146 std::cout
0147 << "PHG4MicromegasHitReco::InitRun\n"
0148 << " m_tmin: " << m_tmin << "ns, m_tmax: " << m_tmax << "ns\n"
0149 << " m_electrons_per_gev: " << m_electrons_per_gev << "\n"
0150 << " m_gain: " << m_gain << "\n"
0151 << " m_cloud_sigma: " << m_cloud_sigma << "cm\n"
0152 << " m_diffusion_trans: " << m_diffusion_trans << "cm/sqrt(cm)\n"
0153 << " m_added_smear_sigma_z: " << m_added_smear_sigma_z << "cm\n"
0154 << " m_added_smear_sigma_rphi: " << m_added_smear_sigma_rphi << "cm\n"
0155 << std::endl;
0156
0157
0158 PHNodeIterator iter(topNode);
0159 auto dstNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "DST"));
0160 if (!dstNode)
0161 {
0162 std::cout << "PHG4MicromegasHitReco::InitRun - DST Node missing, doing nothing." << std::endl;
0163 exit(1);
0164 }
0165
0166
0167 auto hitsetcontainer = findNode::getClass<TrkrHitSetContainer>(topNode, "TRKR_HITSET");
0168 if (!hitsetcontainer)
0169 {
0170 std::cout << "PHG4MicromegasHitReco::InitRun - creating TRKR_HITSET." << std::endl;
0171
0172
0173 PHNodeIterator dstiter(dstNode);
0174 auto trkrnode = dynamic_cast<PHCompositeNode*>(dstiter.findFirst("PHCompositeNode", "TRKR"));
0175 if (!trkrnode)
0176 {
0177 trkrnode = new PHCompositeNode("TRKR");
0178 dstNode->addNode(trkrnode);
0179 }
0180
0181
0182 hitsetcontainer = new TrkrHitSetContainerv1;
0183 auto newNode = new PHIODataNode<PHObject>(hitsetcontainer, "TRKR_HITSET", "PHObject");
0184 trkrnode->addNode(newNode);
0185 }
0186
0187
0188 auto hittruthassoc = findNode::getClass<TrkrHitTruthAssoc>(topNode, "TRKR_HITTRUTHASSOC");
0189 if (!hittruthassoc)
0190 {
0191 std::cout << "PHG4MicromegasHitReco::InitRun - creating TRKR_HITTRUTHASSOC." << std::endl;
0192
0193
0194 PHNodeIterator dstiter(dstNode);
0195 auto trkrnode = dynamic_cast<PHCompositeNode*>(dstiter.findFirst("PHCompositeNode", "TRKR"));
0196 if (!trkrnode)
0197 {
0198 trkrnode = new PHCompositeNode("TRKR");
0199 dstNode->addNode(trkrnode);
0200 }
0201
0202 hittruthassoc = new TrkrHitTruthAssocv1;
0203 auto newNode = new PHIODataNode<PHObject>(hittruthassoc, "TRKR_HITTRUTHASSOC", "PHObject");
0204 trkrnode->addNode(newNode);
0205 }
0206
0207 return Fun4AllReturnCodes::EVENT_OK;
0208 }
0209
0210
0211 int PHG4MicromegasHitReco::process_event(PHCompositeNode* topNode)
0212 {
0213
0214
0215 auto g4hitcontainer = findNode::getClass<PHG4HitContainer>(topNode, "G4HIT_MICROMEGAS");
0216 assert(g4hitcontainer);
0217
0218
0219 m_acts_geometry = findNode::getClass<ActsGeometry>(topNode, "ActsGeometry");
0220 assert(m_acts_geometry);
0221
0222
0223 const auto geonodename = full_geonodename();
0224 auto geonode = findNode::getClass<PHG4CylinderGeomContainer>(topNode, geonodename.c_str());
0225 assert(geonode);
0226
0227
0228 auto trkrhitsetcontainer = findNode::getClass<TrkrHitSetContainer>(topNode, "TRKR_HITSET");
0229 assert(trkrhitsetcontainer);
0230
0231
0232 auto hittruthassoc = findNode::getClass<TrkrHitTruthAssoc>(topNode, "TRKR_HITTRUTHASSOC");
0233 assert(hittruthassoc);
0234
0235
0236 auto layer_range = g4hitcontainer->getLayers();
0237 for (auto layer_it = layer_range.first; layer_it != layer_range.second; ++layer_it)
0238 {
0239
0240 const auto layer = *layer_it;
0241
0242
0243 auto layergeom = dynamic_cast<CylinderGeomMicromegas*>(geonode->GetLayerGeom(layer));
0244 assert(layergeom);
0245
0246
0247
0248
0249
0250
0251
0252
0253 const auto mesh_local_z = layergeom->get_drift_direction() == MicromegasDefs::DriftDirection::OUTWARD ? layergeom->get_thickness() / 2 : -layergeom->get_thickness() / 2;
0254
0255
0256
0257
0258
0259
0260
0261
0262
0263
0264
0265 const PHG4HitContainer::ConstRange g4hit_range = g4hitcontainer->getHits(layer);
0266
0267
0268 for (auto g4hit_it = g4hit_range.first; g4hit_it != g4hit_range.second; ++g4hit_it)
0269 {
0270
0271 PHG4Hit* g4hit = g4hit_it->second;
0272
0273
0274 if (g4hit->get_t(0) > m_tmax)
0275 {
0276 continue;
0277 }
0278 if (g4hit->get_t(1) < m_tmin)
0279 {
0280 continue;
0281 }
0282
0283
0284 TVector3 world_in(g4hit->get_x(0), g4hit->get_y(0), g4hit->get_z(0));
0285 TVector3 world_out(g4hit->get_x(1), g4hit->get_y(1), g4hit->get_z(1));
0286
0287
0288 const int tileid = g4hit->get_property_int(PHG4Hit::prop_index_i);
0289
0290
0291 const auto local_in = layergeom->get_local_from_world_coords(tileid, m_acts_geometry, world_in);
0292 const auto local_out = layergeom->get_local_from_world_coords(tileid, m_acts_geometry, world_out);
0293
0294
0295 const auto nprimary = get_primary_electrons(g4hit);
0296 if (!nprimary)
0297 {
0298 continue;
0299 }
0300
0301
0302 const TrkrDefs::hitsetkey hitsetkey = MicromegasDefs::genHitSetKey(layer, layergeom->get_segmentation_type(), tileid);
0303 const auto hitset_it = trkrhitsetcontainer->findOrAddHitSet(hitsetkey);
0304
0305
0306 using charge_map_t = std::map<int, double>;
0307 charge_map_t total_charges;
0308
0309
0310 for (uint ie = 0; ie < nprimary; ++ie)
0311 {
0312
0313
0314
0315 const auto t = gsl_ran_flat(m_rng.get(), 0.0, 1.0);
0316 auto local = local_in * t + local_out * (1.0 - t);
0317
0318 if (m_diffusion_trans > 0)
0319 {
0320
0321
0322 const double z = local.z();
0323 const double drift_distance = std::abs(z - mesh_local_z);
0324 const double diffusion = gsl_ran_gaussian(m_rng.get(), m_diffusion_trans * std::sqrt(drift_distance));
0325 const double diffusion_angle = gsl_ran_flat(m_rng.get(), -M_PI, M_PI);
0326
0327
0328 local += TVector3(diffusion * std::cos(diffusion_angle), diffusion * std::sin(diffusion_angle), 0);
0329 }
0330
0331 const auto& added_smear_sigma = layergeom->get_segmentation_type() == MicromegasDefs::SegmentationType::SEGMENTATION_PHI ? m_added_smear_sigma_rphi : m_added_smear_sigma_z;
0332
0333 if (added_smear_sigma > 0)
0334 {
0335
0336 const double added_smear_trans = gsl_ran_gaussian(m_rng.get(), added_smear_sigma);
0337 const double added_smear_angle = gsl_ran_flat(m_rng.get(), -M_PI, M_PI);
0338 local += TVector3(added_smear_trans * std::cos(added_smear_angle), added_smear_trans * std::sin(added_smear_angle), 0);
0339 }
0340
0341
0342 const auto fractions = distribute_charge(layergeom, tileid, {local.x(), local.y()}, m_cloud_sigma);
0343
0344
0345 if (Verbosity() > 10)
0346 {
0347 const auto sum = std::accumulate(fractions.begin(), fractions.end(), double(0),
0348 [](double value, const charge_pair_t& pair)
0349 { return value + pair.second; });
0350 std::cout << "PHG4MicromegasHitReco::process_event - sum: " << sum << std::endl;
0351 }
0352
0353
0354 const auto gain = get_single_electron_amplification();
0355
0356
0357 for (const auto& pair : fractions)
0358 {
0359 const int strip = pair.first;
0360 if (strip < 0 || strip >= (int) layergeom->get_strip_count(tileid, m_acts_geometry))
0361 {
0362 continue;
0363 }
0364
0365 const auto it = total_charges.lower_bound(strip);
0366 if (it != total_charges.end() && it->first == strip)
0367 {
0368 it->second += pair.second * gain;
0369 }
0370 else
0371 {
0372 total_charges.insert(it, std::make_pair(strip, pair.second * gain));
0373 }
0374 }
0375 }
0376
0377
0378
0379 for (const auto pair : total_charges)
0380 {
0381
0382 const int strip = pair.first;
0383
0384
0385 TrkrDefs::hitkey hitkey = MicromegasDefs::genHitKey(strip);
0386 auto hit = hitset_it->second->getHit(hitkey);
0387 if (!hit)
0388 {
0389
0390 hit = new TrkrHitv2;
0391 hitset_it->second->addHitSpecificKey(hitkey, hit);
0392 }
0393
0394
0395 hit->addEnergy(pair.second);
0396
0397
0398 hittruthassoc->addAssoc(hitsetkey, hitkey, g4hit_it->first);
0399 }
0400 }
0401 }
0402
0403 return Fun4AllReturnCodes::EVENT_OK;
0404 }
0405
0406
0407 void PHG4MicromegasHitReco::SetDefaultParameters()
0408 {
0409
0410
0411
0412
0413
0414
0415 set_default_double_param("micromegas_tmin", -20);
0416 set_default_double_param("micromegas_tmax", 800);
0417
0418
0419
0420
0421
0422 static constexpr double Ar_dEdx = 2.44;
0423 static constexpr double iC4H10_dEdx = 5.93;
0424 static constexpr double mix_dEdx = 0.9 * Ar_dEdx + 0.1 * iC4H10_dEdx;
0425
0426
0427 static constexpr double Ar_ntot = 94;
0428 static constexpr double iC4H10_ntot = 195;
0429 static constexpr double mix_ntot = 0.9 * Ar_ntot + 0.1 * iC4H10_ntot;
0430
0431
0432 static constexpr double mix_electrons_per_gev = 1e6 * mix_ntot / mix_dEdx;
0433 set_default_double_param("micromegas_electrons_per_gev", mix_electrons_per_gev);
0434
0435
0436 set_default_double_param("micromegas_gain", 2000);
0437
0438
0439 set_default_double_param("micromegas_cloud_sigma", 0.04);
0440
0441
0442 set_default_double_param("micromegas_diffusion_trans", 0.03);
0443
0444
0445 set_default_double_param("micromegas_added_smear_sigma_z", 0);
0446 set_default_double_param("micromegas_added_smear_sigma_rphi", 0);
0447 }
0448
0449
0450 uint PHG4MicromegasHitReco::get_primary_electrons(PHG4Hit* g4hit) const
0451 {
0452 return gsl_ran_poisson(m_rng.get(), g4hit->get_eion() * m_electrons_per_gev);
0453 }
0454
0455
0456 uint PHG4MicromegasHitReco::get_single_electron_amplification() const
0457 {
0458
0459
0460
0461
0462
0463 return gsl_ran_exponential(m_rng.get(), m_gain);
0464 }
0465
0466
0467 PHG4MicromegasHitReco::charge_list_t PHG4MicromegasHitReco::distribute_charge(
0468 CylinderGeomMicromegas* layergeom,
0469 uint tileid,
0470 const TVector2& local_coords,
0471 double sigma) const
0472 {
0473
0474 auto stripnum = layergeom->find_strip_from_local_coords(tileid, m_acts_geometry, local_coords);
0475
0476
0477 if (stripnum < 0)
0478 {
0479 return charge_list_t();
0480 }
0481
0482
0483 const auto pitch = layergeom->get_pitch();
0484
0485
0486 const auto strip_count = layergeom->get_strip_count(tileid, m_acts_geometry);
0487 const int nstrips = 5. * sigma / pitch + 1;
0488 const auto stripnum_min = std::clamp<int>(stripnum - nstrips, 0, strip_count);
0489 const auto stripnum_max = std::clamp<int>(stripnum + nstrips, 0, strip_count);
0490
0491
0492 charge_list_t charge_list;
0493
0494
0495 for (int strip = stripnum_min; strip <= stripnum_max; ++strip)
0496 {
0497
0498 const auto strip_location = layergeom->get_local_coordinates(tileid, m_acts_geometry, strip);
0499
0500
0501
0502
0503
0504
0505 const auto xloc = layergeom->get_segmentation_type() == MicromegasDefs::SegmentationType::SEGMENTATION_PHI ? (strip_location.X() - local_coords.X()) : (strip_location.Y() - local_coords.Y());
0506
0507
0508
0509
0510
0511
0512 const bool zigzag_strips = (layergeom->get_segmentation_type() == MicromegasDefs::SegmentationType::SEGMENTATION_Z);
0513
0514
0515 const auto fraction = zigzag_strips ? get_zigzag_fraction(xloc, sigma, pitch) : get_rectangular_fraction(xloc, sigma, pitch);
0516
0517
0518 charge_list.push_back(std::make_pair(strip, fraction));
0519 }
0520
0521 return charge_list;
0522 }