File indexing completed on 2026-06-05 07:52:43
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053
0054
0055
0056
0057
0058
0059
0060
0061
0062
0063
0064
0065
0066
0067
0068
0069
0070
0071
0072
0073
0074
0075
0076
0077
0078
0079
0080
0081
0082
0083
0084
0085
0086
0087
0088
0089
0090
0091 #include "MPGDTrackerDigi.h"
0092
0093 #include <DD4hep/Alignments.h>
0094 #include <DD4hep/DetElement.h>
0095 #include <DD4hep/IDDescriptor.h>
0096 #include <DD4hep/Objects.h>
0097 #include <DD4hep/Readout.h>
0098 #include <DD4hep/Shapes.h>
0099 #include <DD4hep/VolumeManager.h>
0100 #include <DD4hep/detail/SegmentationsInterna.h>
0101 #include <DDSegmentation/BitFieldCoder.h>
0102 #include <DDSegmentation/CartesianGridUV.h>
0103 #include <DDSegmentation/CartesianGridXY.h>
0104 #include <DDSegmentation/CylindricalGridPhiZ.h>
0105 #include <DDSegmentation/MultiSegmentation.h>
0106 #include <DDSegmentation/Segmentation.h>
0107 #include <Evaluator/DD4hepUnits.h>
0108 #include <Math/GenVector/Cartesian3D.h>
0109 #include <Math/GenVector/DisplacementVector3D.h>
0110 #include <Parsers/Primitives.h>
0111 #include <TGeoMatrix.h>
0112 #include <TMath.h>
0113
0114 #include <algorithms/geo.h>
0115 #include <algorithms/logger.h>
0116 #include <edm4eic/unit_system.h>
0117 #include <edm4hep/EDM4hepVersion.h>
0118 #include <edm4hep/MCParticleCollection.h>
0119 #include <edm4hep/Vector3d.h>
0120 #include <edm4hep/Vector3f.h>
0121 #include <fmt/format.h>
0122 #include <podio/detail/Link.h>
0123 #include <podio/detail/LinkCollectionImpl.h>
0124 #include <algorithm>
0125 #include <array>
0126 #include <cmath>
0127 #include <cstdint>
0128 #include <gsl/pointers>
0129 #include <gsl/util>
0130 #include <initializer_list>
0131 #include <iterator>
0132 #include <map>
0133 #include <memory>
0134 #include <random>
0135 #include <stdexcept>
0136 #include <tuple>
0137 #include <unordered_map>
0138 #include <utility>
0139 #include <vector>
0140
0141 #include "algorithms/digi/MPGDTrackerDigiConfig.h"
0142
0143 using namespace dd4hep;
0144
0145 namespace eicrecon {
0146
0147 void MPGDTrackerDigi::init() {
0148
0149 m_detector = algorithms::GeoSvc::instance().detector();
0150
0151 if (m_cfg.readout.empty()) {
0152 throw std::runtime_error("Readout is empty");
0153 }
0154 try {
0155 m_seg = m_detector->readout(m_cfg.readout).segmentation();
0156 m_id_dec = m_detector->readout(m_cfg.readout).idSpec().decoder();
0157 } catch (const std::runtime_error&) {
0158 critical(R"(Failed to load ID decoder for "{}" readout.)", m_cfg.readout);
0159 throw std::runtime_error("Failed to load ID decoder");
0160 }
0161
0162
0163 parseIDDescriptor();
0164 parseSegmentation();
0165
0166
0167 m_stripRank = [&](CellID vID) {
0168 CellID sID = vID & m_stripBits;
0169 for (int rank = 0; rank < 5; rank++) {
0170 if (sID == m_stripIDs[rank]) {
0171 return rank;
0172 }
0173 }
0174 return -1;
0175 };
0176 m_orientation = [&](CellID vID, CellID vJD) {
0177 int ranki = m_stripRank(vID);
0178 int rankj = m_stripRank(vJD);
0179 if (rankj > ranki) {
0180 return +1;
0181 }
0182 if (rankj < ranki) {
0183 return -1;
0184 }
0185 return 0;
0186 };
0187 m_isUpstream = [](int orientation, unsigned int status) {
0188
0189 bool isUpstream =
0190 (orientation < 0 && (status & 0x2)) ||
0191 (orientation > 0 && (status & 0x8)) ||
0192 (orientation == 0 && (status & 0x102) == 0x102);
0193 return isUpstream;
0194 };
0195 m_isDownstream = [](int orientation, unsigned int status) {
0196
0197 bool isDownstream =
0198 (orientation > 0 && (status & 0x1)) ||
0199 (orientation < 0 && (status & 0x4)) ||
0200 (orientation == 0 && (status & 0x101) == 0x101);
0201 return isDownstream;
0202 };
0203
0204 m_toleranceFactor = [](double P) {
0205 int factor = 0;
0206 if (P < 1 * dd4hep::MeV) {
0207 factor = 4;
0208 } else if (P < 10 * dd4hep::MeV) {
0209 factor = 2;
0210 } else {
0211 factor = 1;
0212 }
0213 return factor;
0214 };
0215
0216
0217 m_binToPosition = [](FieldID bin, double cellSize, double offset) {
0218 return bin * cellSize + offset;
0219 };
0220 }
0221
0222
0223 void getLocalPosMom(const edm4hep::SimTrackerHit& sim_hit, const TGeoHMatrix& toModule,
0224 double* lpos, double* lmom);
0225 bool cExtrapolate(const double* lpos, const double* lmom,
0226 double rT,
0227 double* lext);
0228 double getRef2Cur(DetElement refVol, DetElement curVol);
0229 bool bExtrapolate(const double* lpos, const double* lmom,
0230 double zT,
0231 double* lext);
0232 std::string inconsistency(const edm4hep::EventHeader& event, unsigned int status, CellID cID,
0233 const double* lpos, const double* lmom);
0234 std::string oddity(const edm4hep::EventHeader& event, unsigned int status, double dist, CellID cID,
0235 const double* lpos, const double* lmom, CellID cJD, const double* lpoj,
0236 const double* lmoj);
0237 double outInDistance(int shape, int orientation, double lintos[][3], double louts[][3],
0238 double* lmom, double* lmoj);
0239 void flagUnexpected(const edm4hep::EventHeader& event, int shape, double expected,
0240 const edm4hep::SimTrackerHit& sim_hit, double* lpini, double* lpend,
0241 double* lpos, double* lmom);
0242
0243 void MPGDTrackerDigi::process(const MPGDTrackerDigi::Input& input,
0244 const MPGDTrackerDigi::Output& output) const {
0245
0246 const auto [headers, sim_hits] = input;
0247 #if EDM4EIC_BUILD_VERSION >= EDM4EIC_VERSION(8, 7, 0)
0248 auto [raw_hits, links, associations] = output;
0249 #else
0250 auto [raw_hits, associations] = output;
0251 #endif
0252
0253
0254 auto seed = m_uid.getUniqueID(*headers, name());
0255 std::default_random_engine generator(seed);
0256 std::normal_distribution<double> gaussian;
0257
0258
0259 std::unordered_map<std::uint64_t, edm4eic::MutableRawTrackerHit> cell_hit_maps[2];
0260
0261 std::map<std::uint64_t, std::vector<std::uint64_t>> stripID2cIDs;
0262
0263 const VolumeManager& volman = m_detector->volumeManager();
0264
0265
0266
0267
0268 const edm4hep::EventHeader& header = headers->at(0);
0269
0270
0271 size_t sim_size = sim_hits->size();
0272 for (int idx = 0; idx < (int)sim_size; idx++) {
0273 const edm4hep::SimTrackerHit& sim_hit = sim_hits->at(idx);
0274
0275
0276
0277
0278
0279
0280
0281 double time_smearing = gaussian(generator) * m_cfg.timeResolution;
0282
0283
0284 CellID refID = sim_hit.getCellID() & m_moduleBits;
0285 DetElement refVol = volman.lookupDetElement(refID);
0286
0287
0288
0289
0290
0291 double lpos[3], eDep, time;
0292 std::vector<std::uint64_t> cIDs;
0293 const auto& shape = refVol.solid();
0294 if (std::string_view{shape.type()} == "TGeoTubeSeg") {
0295
0296 if (!cCoalesceExtend(input, idx, cIDs, lpos, eDep, time))
0297 continue;
0298 } else if (std::string_view{shape.type()} == "TGeoBBox") {
0299
0300 if (!bCoalesceExtend(input, idx, cIDs, lpos, eDep, time))
0301 continue;
0302 } else {
0303 critical(R"(Bad input data: CellID {:x} has invalid shape "{}")", refID, shape.type());
0304 throw std::runtime_error(R"(Inconsistency: Inappropriate SimHits fed to "MPGDTrackerDigi".)");
0305 }
0306
0307
0308 double surfPos[2];
0309 Position locPos;
0310 if (std::string_view{shape.type()} == "TGeoTubeSeg") {
0311
0312 const Tube& tRef = refVol.solid();
0313 double R = (tRef.rMin() + tRef.rMax()) / 2;
0314 double phi = atan2(lpos[1], lpos[0]);
0315 surfPos[0] = phi * R;
0316 surfPos[1] = lpos[2];
0317 locPos = Position(R * cos(phi), R * sin(phi), lpos[2]);
0318 } else {
0319 locPos = Position(lpos[0], lpos[1], lpos[2]);
0320 if (m_gridAngle != 0.0) {
0321 dd4hep::DDSegmentation::Vector3D position;
0322 position.X = lpos[0];
0323 position.Y = lpos[1];
0324 position = RotationZ(m_gridAngle) * position;
0325 surfPos[0] = position.X;
0326 surfPos[1] = position.Y;
0327 } else {
0328 surfPos[0] = lpos[0];
0329 surfPos[1] = lpos[1];
0330 }
0331 }
0332
0333
0334 for (int pn = 0; pn < 2; pn++) {
0335
0336
0337 Cluster cluster;
0338 int status = get2HitCluster(refID, locPos, surfPos, pn, generator, cluster);
0339 if (status != 0) {
0340 CellID vIDs = (std::uint32_t)refID;
0341 CellID hIDs = refID >> 32;
0342 error(R"(SimHit (= 0x{:08x}, 0x{:08x}, {:.2f},{:.2f} cm) beyond limits of {}Strips.)", hIDs,
0343 vIDs, surfPos[0] / cm, surfPos[1] / cm, (pn != 0) ? 'n' : 'p');
0344 }
0345
0346 if (level() >= algorithms::LogLevel::kDebug) {
0347 std::string sCellID = (pn != 0) ? "cellIDn" : "cellIDp";
0348 if (pn == 0) {
0349 debug(" =>=>");
0350 }
0351 for (auto clusterHit : cluster) {
0352 CellID cID = clusterHit.first;
0353 CellID hID = cID >> 32;
0354 CellID vID = (std::uint32_t)cID;
0355 debug("Hit {} = 0x{:08x}, 0x{:08x}", sCellID, hID, vID);
0356 debug(" edep = {:.0f} [eV]", clusterHit.second * eDep / eV);
0357 }
0358 }
0359
0360
0361
0362 std::unordered_map<std::uint64_t, edm4eic::MutableRawTrackerHit>& cell_hit_map =
0363 gsl::at(cell_hit_maps, pn);
0364 double g = m_cfg.gain;
0365 for (auto clusterHit : cluster) {
0366 CellID cID = clusterHit.first;
0367 double f = clusterHit.second;
0368
0369
0370
0371
0372
0373 if (eDep < m_cfg.threshold) {
0374 debug(" eDep {:.2f} is below threshold of {:.2f} [keV]", eDep, m_cfg.threshold / keV);
0375 continue;
0376 }
0377 stripID2cIDs[cID] = cIDs;
0378 double result_time = time + time_smearing;
0379 auto hit_time_stamp = (std::int32_t)(result_time * 1e3);
0380 if (!cell_hit_map.contains(cID)) {
0381
0382 cell_hit_map[cID] = {
0383 cID, (std::int32_t)std::llround(eDep * 1e6 * f * g),
0384 hit_time_stamp
0385 };
0386 } else {
0387
0388 auto& hit = cell_hit_map[cID];
0389 debug(" Hit already exists in cell ID={}, prev. hit time: {}", cID, hit.getTimeStamp());
0390
0391
0392 hit.setTimeStamp(std::min(hit_time_stamp, hit.getTimeStamp()));
0393
0394
0395 auto charge = hit.getCharge();
0396 hit.setCharge(charge + (std::int32_t)std::llround(eDep * 1e6 * f * g));
0397 }
0398 }
0399 }
0400 }
0401
0402
0403 for (auto& cell_hit_map : cell_hit_maps) {
0404 for (auto item : cell_hit_map) {
0405 raw_hits->push_back(item.second);
0406 CellID stripID = item.first;
0407 const auto is = stripID2cIDs.find(stripID);
0408 if (is == stripID2cIDs.end()) {
0409 error(R"(Inconsistency: CellID {:x} not found in "stripID2cIDs" map)", stripID);
0410 throw std::runtime_error(R"(Inconsistency in the handling of "stripID2cIDs" map)");
0411 }
0412 std::vector<std::uint64_t> cIDs = is->second;
0413 for (CellID cID : cIDs) {
0414 for (const auto& sim_hit : *sim_hits) {
0415 if (sim_hit.getCellID() == cID) {
0416 #if EDM4EIC_BUILD_VERSION >= EDM4EIC_VERSION(8, 7, 0)
0417
0418 auto link = links->create();
0419 link.setFrom(item.second);
0420 link.setTo(sim_hit);
0421 link.setWeight(1.0);
0422 #endif
0423
0424 auto hitassoc = associations->create();
0425 hitassoc.setWeight(1.0);
0426 hitassoc.setRawHit(item.second);
0427 hitassoc.setSimHit(sim_hit);
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 void MPGDTrackerDigi::parseIDDescriptor() {
0463
0464
0465
0466
0467
0468
0469
0470 debug(R"(Parsing IDDescriptor for "{}" readout)", m_cfg.readout);
0471 CellID sensorBits = 0;
0472 for (int field = 0; field < 5; field++) {
0473 const char* fieldName = m_fieldNames[field];
0474 CellID fieldID = 0;
0475 try {
0476 fieldID = m_id_dec->get(~((CellID)0x0), fieldName);
0477 } catch (const std::runtime_error& error) {
0478 critical(R"(No field "{}" in IDDescriptor of readout "{}".)", fieldName, m_cfg.readout);
0479 throw std::runtime_error("Invalid IDDescriptor");
0480 }
0481 const BitFieldElement& fieldElement = (*m_id_dec)[fieldName];
0482 CellID fieldBits = fieldID << fieldElement.offset();
0483 m_volumeBits |= fieldBits;
0484
0485 if (std::string_view{fieldName} == "strip") {
0486 m_stripBits = fieldBits;
0487
0488 CellID bits[5] = {0x3, 0x1, 0, 0x2, 0x4};
0489 for (int subVolume = 0; subVolume < 5; subVolume++) {
0490 m_stripIDs[subVolume] = bits[subVolume] << fieldElement.offset();
0491 }
0492 }
0493
0494 else if (m_cfg.readout == "MPGDBarrelHits" && std::string_view{fieldName} == "sensor") {
0495 sensorBits = fieldBits;
0496 m_sensorOffset = fieldElement.offset();
0497 }
0498 }
0499
0500 m_moduleBits = m_volumeBits & ~m_stripBits;
0501 m_pStripBit = m_stripIDs[1];
0502 m_nStripBit = m_stripIDs[3];
0503
0504 if (m_cfg.readout == "MPGDBarrelHits") {
0505 m_sensorStripBits = sensorBits | m_stripBits;
0506 } else {
0507 m_sensorStripBits = m_stripBits;
0508 }
0509
0510
0511
0512
0513 int coordOffsets[2] = {64, 64};
0514 for (int pn = 0; pn < 2; pn++) {
0515 std::string coordName;
0516 if (m_cfg.readout == "MPGDBarrelHits") {
0517 coordName = pn ? "z" : "phi";
0518 } else if (m_cfg.readout == "OuterMPGDBarrelHits") {
0519 coordName = pn ? "v" : "u";
0520 } else {
0521 coordName = pn ? "y" : "x";
0522 }
0523 try {
0524 m_stripIndices[pn] = m_id_dec->index(coordName);
0525 } catch (const std::runtime_error& error) {
0526 critical(R"(No "{}" field in IDDescriptor of readout "{}".)", coordName, m_cfg.readout);
0527 throw std::runtime_error("Invalid IDDescriptor");
0528 }
0529 const BitFieldElement& fieldElement = (*m_id_dec)[coordName];
0530 int offset = fieldElement.offset();
0531 m_stripIncs[pn] = ((CellID)0x1) << offset;
0532 if (offset < coordOffsets[pn])
0533 coordOffsets[pn] = offset;
0534 }
0535 if (coordOffsets[0] != 32 || coordOffsets[1] != 48) {
0536 critical(R"(Coordinate fields in IDDescriptor of readout "{}" do not start @ bits 32 and 48.)",
0537 m_cfg.readout);
0538 throw std::runtime_error("Invalid IDDescriptor");
0539 }
0540 }
0541 void MPGDTrackerDigi::parseSegmentation() {
0542 debug(R"(Find valid "MultiSegmentation" for "{}" readout.)", m_cfg.readout);
0543
0544
0545 std::function<void(int, double)> checkResolutionVsPitch = [&](int pn, double pitch) {
0546 double sigma = m_cfg.stripResolutions[pn];
0547 if (m_truncation * sigma > pitch) {
0548 critical(R"(stripResolutions[{}] (= {}) too large for pitch (={}) of "{}" readout.)", pn,
0549 sigma, pitch, m_cfg.readout);
0550 throw std::runtime_error("Space resolution parameter too large");
0551 }
0552 return;
0553 };
0554 using Segmentation = dd4hep::DDSegmentation::Segmentation;
0555 const Segmentation* segmentation = m_seg->segmentation;
0556
0557 unsigned int required = 0, fulfilled = 0;
0558 if (segmentation->type() == "MultiSegmentation") {
0559 using MultiSegmentation = dd4hep::DDSegmentation::MultiSegmentation;
0560 const auto* multiSeg = dynamic_cast<const MultiSegmentation*>(segmentation);
0561 StripParameters pars;
0562 if (m_cfg.readout == "MPGDBarrelHits") {
0563
0564 required = 0x3 | m_pStripBit | m_nStripBit;
0565 unsigned int innerBit = 0x0, outerBit = 0x1;
0566 for (unsigned int sectorBit : {innerBit, outerBit}) {
0567 CellID sensorID = ((CellID)sectorBit) << m_sensorOffset;
0568 const Segmentation& sectorSeg = multiSeg->subsegmentation(sensorID);
0569 if (multiSeg->type() != "MultiSegmentation")
0570 continue;
0571 const auto& stripMultiSeg = dynamic_cast<const MultiSegmentation&>(sectorSeg);
0572 for (CellID stripID : {m_pStripBit, m_nStripBit}) {
0573 const Segmentation& stripSeg = stripMultiSeg.subsegmentation(stripID);
0574 if (stripSeg.type() != "CylindricalGridPhiZ") {
0575 critical(
0576 R"(Segmentation type for "{}" readout = "{}", whereas expected = "CylindricalGridPhiZ".)",
0577 m_cfg.readout, stripSeg.type());
0578 continue;
0579 }
0580 const dd4hep::DDSegmentation::CylindricalGridPhiZ& gridPhiZ =
0581 dynamic_cast<const dd4hep::DDSegmentation::CylindricalGridPhiZ&>(stripSeg);
0582 fulfilled |= sectorBit == innerBit ? 0x1 : 0x2;
0583 fulfilled |= stripID;
0584
0585 double radius = gridPhiZ.radius();
0586 int pn;
0587 std::string stripName;
0588 if (stripID == m_pStripBit) {
0589 pars.pitch = gridPhiZ.gridSizePhi() * radius;
0590 pars.offset = gridPhiZ.offsetPhi() * radius;
0591 pn = 0;
0592 } else {
0593 pars.pitch = gridPhiZ.gridSizeZ();
0594 pars.offset = gridPhiZ.offsetZ();
0595 pn = 1;
0596 }
0597 pars.index = m_stripIndices[pn];
0598
0599 checkResolutionVsPitch(pn, pars.pitch);
0600 int nStrips = m_cfg.stripNumbers[pn];
0601
0602
0603
0604 pars.max = (nStrips / 2 - .5) * pars.pitch + pars.offset;
0605 pars.min = (-nStrips / 2 - .5) * pars.pitch + pars.offset;
0606
0607 CellID sensorStripID = sensorID | stripID;
0608 m_stripParameters[sensorStripID] = pars;
0609 }
0610 }
0611 } else if (m_cfg.readout == "OuterMPGDBarrelHits") {
0612
0613 required = m_pStripBit | m_nStripBit;
0614 double gridAngles[2];
0615 for (unsigned int stripID : {m_pStripBit, m_nStripBit}) {
0616 const dd4hep::DDSegmentation::Segmentation& stripSeg = multiSeg->subsegmentation(stripID);
0617 if (stripSeg.type() != "CartesianGridUV") {
0618 critical(
0619 R"(Segmentation type for "{}" readout = "{}", whereas expected = "CartesianGridUV".)",
0620 m_cfg.readout, stripSeg.type());
0621 continue;
0622 }
0623 const dd4hep::DDSegmentation::CartesianGridUV& gridUV =
0624 dynamic_cast<const dd4hep::DDSegmentation::CartesianGridUV&>(stripSeg);
0625 fulfilled |= stripID;
0626
0627 int pn;
0628 std::string stripName;
0629 if (stripID == m_pStripBit) {
0630 pars.pitch = gridUV.gridSizeU();
0631 pars.offset = gridUV.offsetU();
0632 pn = 0;
0633 } else {
0634 pars.pitch = gridUV.gridSizeV();
0635 pars.offset = gridUV.offsetV();
0636 pn = 1;
0637 }
0638 pars.index = m_stripIndices[pn];
0639 gridAngles[pn] = gridUV.gridAngle();
0640
0641 checkResolutionVsPitch(pn, pars.pitch);
0642 int nStrips = m_cfg.stripNumbers[pn];
0643 pars.max = (nStrips / 2 - .5) * pars.pitch + pars.offset;
0644 pars.min = (-nStrips / 2 - .5) * pars.pitch + pars.offset;
0645
0646 m_stripParameters[stripID] = pars;
0647 }
0648 if (fabs(gridAngles[0] - gridAngles[1]) > 1e-6) {
0649 critical(
0650 R"(Inconsistent gridAngles of "CartesianGridUV" for readout "{}": 'p' = {:.6f}, 'n' = {:.6f}.)",
0651 m_cfg.readout, gridAngles[0], gridAngles[1]);
0652 fulfilled = 0;
0653 }
0654 m_gridAngle = gridAngles[0];
0655 } else {
0656
0657 required = m_pStripBit | m_nStripBit;
0658 for (unsigned int stripID : {m_pStripBit, m_nStripBit}) {
0659 const dd4hep::DDSegmentation::Segmentation& stripSeg = multiSeg->subsegmentation(stripID);
0660 if (stripSeg.type() != "CartesianGridXY") {
0661 critical(
0662 R"(Segmentation type for "{}" readout = "{}", whereas expected = "CartesianGridXY".)",
0663 m_cfg.readout, stripSeg.type());
0664 continue;
0665 }
0666 const dd4hep::DDSegmentation::CartesianGridXY& gridXY =
0667 dynamic_cast<const dd4hep::DDSegmentation::CartesianGridXY&>(stripSeg);
0668 fulfilled |= stripID;
0669
0670 int pn;
0671 std::string stripName;
0672 if (stripID == m_pStripBit) {
0673 pars.pitch = gridXY.gridSizeX();
0674 pars.offset = gridXY.offsetX();
0675 pn = 0;
0676 } else {
0677 pars.pitch = gridXY.gridSizeY();
0678 pars.offset = gridXY.offsetY();
0679 pn = 1;
0680 }
0681 pars.index = m_stripIndices[pn];
0682
0683 checkResolutionVsPitch(pn, pars.pitch);
0684 int nStrips = m_cfg.stripNumbers[pn];
0685 pars.max = (nStrips / 2 - .5) * pars.pitch + pars.offset;
0686 pars.min = (-nStrips / 2 - .5) * pars.pitch + pars.offset;
0687
0688 m_stripParameters[stripID] = pars;
0689 }
0690 }
0691 } else {
0692 critical(
0693 R"(Segmentation type for "{}" readout = "{}", whereas expected = "MultiSegmentation".)",
0694 m_cfg.readout, segmentation->type());
0695 }
0696 if (fulfilled != required) {
0697 critical(R"(Error retrieving Segmentation parameters for "{}" readout.)", m_cfg.readout);
0698 throw std::runtime_error("Error retrieving Segmentation parameters");
0699 }
0700 }
0701
0702
0703
0704
0705
0706
0707
0708
0709
0710
0711
0712
0713
0714
0715 bool MPGDTrackerDigi::cCoalesceExtend(const Input& input, int& idx,
0716 std::vector<std::uint64_t>& cIDs, double* lpos, double& eDep,
0717 double& time) const {
0718 const auto [headers, sim_hits] = input;
0719 const edm4hep::EventHeader& header = headers->at(0);
0720 const edm4hep::SimTrackerHit& sim_hit = sim_hits->at(idx);
0721 CellID vID = sim_hit.getCellID() & m_volumeBits;
0722 CellID refID = vID & m_moduleBits;
0723 const VolumeManager& volman = m_detector->volumeManager();
0724 DetElement refVol = volman.lookupDetElement(refID);
0725
0726
0727 const TGeoHMatrix toRefVol = refVol.nominal().worldTransformation();
0728 double lmom[3];
0729 getLocalPosMom(sim_hit, toRefVol, lpos, lmom);
0730 const double edmm = edm4eic::unit::mm, ed2dd = dd4hep::mm / edmm;
0731 using dd4hep::mm;
0732 using edm4eic::unit::eV, edm4eic::unit::GeV;
0733
0734 eDep = sim_hit.getEDep();
0735 time = sim_hit.getTime();
0736
0737 const Tube& tRef = refVol.solid();
0738 double dZ = tRef.dZ();
0739
0740
0741
0742
0743
0744 double startPhi = tRef.startPhi() * radian;
0745 startPhi -= 2 * TMath::Pi();
0746 double endPhi = tRef.endPhi() * radian;
0747 endPhi -= 2 * TMath::Pi();
0748
0749 DetElement curVol = volman.lookupDetElement(vID);
0750 const Tube& tCur = curVol.solid();
0751 double rMin = tCur.rMin(), rMax = tCur.rMax();
0752
0753 double lintos[2][3], louts[2][3], lpini[3], lpend[3], lmend[3];
0754 std::copy(std::begin(lmom), std::end(lmom), std::begin(lmend));
0755 unsigned int status =
0756 cTraversing(lpos, lmom, sim_hit.getPathLength() * ed2dd, sim_hit.isProducedBySecondary(),
0757 rMin, rMax, dZ, startPhi, endPhi, lintos, louts, lpini, lpend);
0758 if (status & m_inconsistency) {
0759 error(inconsistency(header, status, sim_hit.getCellID(), lpos, lmom));
0760 return false;
0761 }
0762 cIDs.push_back(sim_hit.getCellID());
0763 std::vector<int> subHitList;
0764 if (level() >= algorithms::LogLevel::kDebug) {
0765 subHitList.push_back(idx);
0766 }
0767
0768 bool isContinuation = status & (m_intoLower | m_intoUpper);
0769 bool hasContinuation = status & (m_outLower | m_outUpper);
0770 bool canReEnter = status & m_canReEnter;
0771 int rank = m_stripRank(vID);
0772 if (!canReEnter) {
0773 if (rank == 0 && (status & m_intoLower))
0774 isContinuation = false;
0775 if (rank == 0 && (status & m_outLower))
0776 hasContinuation = false;
0777 }
0778 if (rank == 4 && (status & m_intoUpper))
0779 isContinuation = false;
0780 if (rank == 4 && (status & m_outUpper))
0781 hasContinuation = false;
0782 if (hasContinuation) {
0783
0784 int jdx;
0785 CellID vIDPrv = vID;
0786 size_t sim_size = sim_hits->size();
0787 for (jdx = idx + 1; jdx < (int)sim_size; jdx++) {
0788 const edm4hep::SimTrackerHit& sim_hjt = sim_hits->at(jdx);
0789 CellID vJD = sim_hjt.getCellID() & m_volumeBits;
0790
0791
0792 int orientation = m_orientation(vIDPrv, vJD);
0793 bool isUpstream = m_isUpstream(orientation, status);
0794 bool pmoStatus = samePMO(sim_hit, sim_hjt);
0795 if (!pmoStatus || !isUpstream) {
0796 if ((pmoStatus && !isUpstream) && !sim_hit.isProducedBySecondary()) {
0797
0798
0799 double P = sqrt(lmom[0] * lmom[0] + lmom[1] * lmom[1] + lmom[2] * lmom[2]);
0800 if (P > 10 * dd4hep::MeV)
0801 debug(inconsistency(header, 0, sim_hit.getCellID(), lpos, lmom));
0802 }
0803 break;
0804 }
0805
0806 curVol = volman.lookupDetElement(vJD);
0807 const Tube& tubj = curVol.solid();
0808 rMin = tubj.rMin();
0809 rMax = tubj.rMax();
0810 double lpoj[3], lmoj[3];
0811 getLocalPosMom(sim_hjt, toRefVol, lpoj, lmoj);
0812
0813 double ljns[2][3], lovts[2][3], lpjni[3], lpfnd[3];
0814 status =
0815 cTraversing(lpoj, lmoj, sim_hjt.getPathLength() * ed2dd, sim_hit.isProducedBySecondary(),
0816 rMin, rMax, dZ, startPhi, endPhi, ljns, lovts, lpjni, lpfnd);
0817 if (status & m_inconsistency) {
0818 error(inconsistency(header, status, sim_hjt.getCellID(), lpoj, lmoj));
0819 break;
0820 }
0821
0822 bool jsDownstream = m_isDownstream(orientation, status);
0823 if (!jsDownstream)
0824 break;
0825
0826 double dist = outInDistance(0, orientation, ljns, louts, lmom, lmoj);
0827
0828 double P = sqrt(lmom[0] * lmom[0] + lmom[1] * lmom[1] + lmom[2] * lmom[2]);
0829 double tolerance = m_toleranceFactor(P) * 25 * dd4hep::um;
0830 bool isCompatible = dist > 0 && dist < tolerance;
0831 if (!isCompatible) {
0832 if (!sim_hit.isProducedBySecondary())
0833 debug(oddity(header, status, dist, sim_hit.getCellID(), lpos, lmom,
0834 sim_hjt.getCellID(), lpoj, lmoj));
0835 break;
0836 }
0837
0838 vIDPrv = vJD;
0839 eDep += sim_hjt.getEDep();
0840 for (int i = 0; i < 3; i++) {
0841 lpend[i] = lpfnd[i];
0842 lmend[i] = lmoj[i];
0843 }
0844
0845 cIDs.push_back(sim_hjt.getCellID());
0846 if (level() >= algorithms::LogLevel::kDebug) {
0847 subHitList.push_back(jdx);
0848 }
0849
0850 hasContinuation = status & 0xa;
0851 canReEnter = status & 0x100;
0852 if (!canReEnter && m_stripRank(vJD) == 4)
0853 hasContinuation = false;
0854 if (!hasContinuation) {
0855 jdx++;
0856 break;
0857 } else {
0858 for (int i = 0; i < 3; i++) {
0859 louts[0][i] = lovts[0][i];
0860 louts[1][i] = lovts[1][i];
0861 }
0862 }
0863 }
0864 idx = jdx - 1;
0865 }
0866
0867 if (sim_hit.isProducedBySecondary() && cIDs.size() < 2)
0868 if (denyExtension(sim_hit, tCur.rMax() - tCur.rMin())) {
0869 isContinuation = hasContinuation = false;
0870 }
0871 for (int io = 0; io < 2; io++) {
0872 if ((io == 0 && !isContinuation) || (io == 1 && !hasContinuation))
0873 continue;
0874 int direction = io ? +1 : -1;
0875 extendHit(refID, cIDs, direction, lpini, lmom, lpend, lmend);
0876 }
0877
0878 flagUnexpected(header, 0, (tRef.rMin() + tRef.rMax()) / 2, sim_hit, lpini, lpend, lpos, lmom);
0879
0880 double DoF2 = 0, dir = 0;
0881 for (int i = 0; i < 3; i++) {
0882 double neu = (lpini[i] + lpend[i]) / 2, alt = lpos[i];
0883 lpos[i] = neu;
0884 double d = neu - alt;
0885 dir += d * lmom[i];
0886 DoF2 += d * d;
0887 }
0888
0889 time += ((dir > 0) ? 1 : ((dir < 0) ? -1 : 0)) * sqrt(DoF2) / dd4hep::c_light;
0890 if (level() >= algorithms::LogLevel::kDebug) {
0891 debug("--------------------");
0892 printSubHitList(input, subHitList);
0893 debug(" =");
0894
0895 Position locPos(lpos[0], lpos[1], lpos[2]);
0896 Position globPos = refVol.nominal().localToWorld(locPos);
0897 debug(" position = ({:7.2f},{:7.2f},{:7.2f}) [mm]", globPos.X() / mm, globPos.Y() / mm,
0898 globPos.Z() / mm);
0899 debug(" edep = {:.0f} [eV]", eDep / eV);
0900 debug(" time = {:.2f} [ns]", time);
0901 }
0902 return true;
0903 }
0904 bool MPGDTrackerDigi::bCoalesceExtend(const Input& input, int& idx,
0905 std::vector<std::uint64_t>& cIDs, double* lpos, double& eDep,
0906 double& time) const {
0907 const auto [headers, sim_hits] = input;
0908 const edm4hep::EventHeader& header = headers->at(0);
0909 const edm4hep::SimTrackerHit& sim_hit = sim_hits->at(idx);
0910 CellID vID = sim_hit.getCellID() & m_volumeBits;
0911 CellID refID = vID & m_moduleBits;
0912 const VolumeManager& volman = m_detector->volumeManager();
0913 DetElement refVol = volman.lookupDetElement(refID);
0914
0915
0916 const TGeoHMatrix toRefVol = refVol.nominal().worldTransformation();
0917 double lmom[3];
0918 getLocalPosMom(sim_hit, toRefVol, lpos, lmom);
0919 const double edmm = edm4eic::unit::mm, ed2dd = dd4hep::mm / edmm;
0920 using dd4hep::mm;
0921 using edm4eic::unit::eV, edm4eic::unit::GeV;
0922
0923 eDep = sim_hit.getEDep();
0924 time = sim_hit.getTime();
0925
0926 const Box& bRef = refVol.solid();
0927 double dX = bRef.x(), dY = bRef.y();
0928
0929 DetElement curVol = volman.lookupDetElement(vID);
0930 const Box& bCur = curVol.solid();
0931 double dZ = bCur.z();
0932 double ref2Cur = getRef2Cur(refVol, curVol);
0933
0934 double lintos[2][3], louts[2][3], lpini[3], lpend[3], lmend[3];
0935 std::copy(std::begin(lmom), std::end(lmom), std::begin(lmend));
0936 unsigned int status =
0937 bTraversing(lpos, lmom, ref2Cur, sim_hit.getPathLength() * ed2dd,
0938 sim_hit.isProducedBySecondary(), dZ, dX, dY, lintos, louts, lpini, lpend);
0939 if (status & m_inconsistency) {
0940 error(inconsistency(header, status, sim_hit.getCellID(), lpos, lmom));
0941 return false;
0942 }
0943 cIDs.push_back(sim_hit.getCellID());
0944 std::vector<int> subHitList;
0945 if (level() >= algorithms::LogLevel::kDebug) {
0946 subHitList.push_back(idx);
0947 }
0948
0949 int rank = m_stripRank(vID);
0950 bool isContinuation = status & (m_intoLower | m_intoUpper);
0951 bool hasContinuation = status & (m_outLower | m_outUpper);
0952 if ((rank == 0 && (status & m_intoLower)) || (rank == 4 && (status & m_intoUpper)))
0953 isContinuation = false;
0954 if ((rank == 0 && (status & m_outLower)) || (rank == 4 && (status & m_outUpper)))
0955 hasContinuation = false;
0956 if (hasContinuation) {
0957
0958 int jdx;
0959 CellID vIDPrv = vID;
0960 size_t sim_size = sim_hits->size();
0961 for (jdx = idx + 1; jdx < (int)sim_size; jdx++) {
0962 const edm4hep::SimTrackerHit& sim_hjt = sim_hits->at(jdx);
0963 CellID vJD = sim_hjt.getCellID() & m_volumeBits;
0964 int orientation = m_orientation(vIDPrv, vJD);
0965 bool isUpstream = m_isUpstream(orientation, status);
0966 bool pmoStatus = samePMO(sim_hit, sim_hjt);
0967 if (!pmoStatus || !isUpstream) {
0968 if ((pmoStatus && !isUpstream) && !sim_hit.isProducedBySecondary()) {
0969
0970 double P = sqrt(lmom[0] * lmom[0] + lmom[1] * lmom[1] + lmom[2] * lmom[2]);
0971 if (P > 10 * dd4hep::MeV)
0972 debug(inconsistency(header, 0, sim_hit.getCellID(), lpos, lmom));
0973 }
0974 break;
0975 }
0976
0977 curVol = volman.lookupDetElement(vJD);
0978 const Box& boxj = curVol.solid();
0979 dZ = boxj.z();
0980 double ref2j = getRef2Cur(refVol, curVol);
0981
0982 double lpoj[3], lmoj[3];
0983 getLocalPosMom(sim_hjt, toRefVol, lpoj, lmoj);
0984 double ljns[2][3], lovts[2][3], lpjni[3], lpfnd[3];
0985 status = bTraversing(lpoj, lmoj, ref2j, sim_hjt.getPathLength() * ed2dd,
0986 sim_hit.isProducedBySecondary(), dZ, dX, dY, ljns, lovts, lpjni, lpfnd);
0987 if (status & m_inconsistency) {
0988 error(inconsistency(header, status, sim_hjt.getCellID(), lpoj, lmoj));
0989 break;
0990 }
0991
0992 bool jsDownstream = m_isDownstream(orientation, status);
0993 if (!jsDownstream)
0994 break;
0995
0996 double dist = outInDistance(1, orientation, ljns, louts, lmom, lmoj);
0997
0998 double P = sqrt(lmom[0] * lmom[0] + lmom[1] * lmom[1] + lmom[2] * lmom[2]);
0999 double tolerance = m_toleranceFactor(P) * 25 * dd4hep::um;
1000 bool isCompatible = dist > 0 && dist < tolerance;
1001 if (!isCompatible) {
1002 if (!sim_hit.isProducedBySecondary())
1003 debug(oddity(header, status, dist, sim_hit.getCellID(), lpos, lmom,
1004 sim_hjt.getCellID(), lpoj, lmoj));
1005 break;
1006 }
1007
1008 vIDPrv = vJD;
1009 eDep += sim_hjt.getEDep();
1010 for (int i = 0; i < 3; i++) {
1011 lpend[i] = lpfnd[i];
1012 lmend[i] = lmoj[i];
1013 }
1014
1015 cIDs.push_back(sim_hjt.getCellID());
1016 if (level() >= algorithms::LogLevel::kDebug) {
1017 subHitList.push_back(jdx);
1018 }
1019
1020 hasContinuation = status & 0xa;
1021 if (!hasContinuation) {
1022 jdx++;
1023 break;
1024 } else {
1025 for (int i = 0; i < 3; i++) {
1026 louts[0][i] = lovts[0][i];
1027 louts[1][i] = lovts[1][i];
1028 }
1029 }
1030 }
1031 idx = jdx - 1;
1032 }
1033
1034 if (sim_hit.isProducedBySecondary() && cIDs.size() < 2)
1035 if (denyExtension(sim_hit, bCur.z())) {
1036 isContinuation = hasContinuation = false;
1037 }
1038 for (int io = 0; io < 2; io++) {
1039 if ((io == 0 && !isContinuation) || (io == 1 && !hasContinuation))
1040 continue;
1041 int direction = io ? +1 : -1;
1042 extendHit(refID, cIDs, direction, lpini, lmom, lpend, lmend);
1043 }
1044
1045 flagUnexpected(header, 1, 0, sim_hit, lpini, lpend, lpos, lmom);
1046
1047 double DoF2 = 0, dir = 0;
1048 for (int i = 0; i < 3; i++) {
1049 double neu = (lpini[i] + lpend[i]) / 2, alt = lpos[i];
1050 lpos[i] = neu;
1051 double d = neu - alt;
1052 dir += d * lmom[i];
1053 DoF2 += d * d;
1054 }
1055
1056 time += ((dir > 0) ? 1 : ((dir < 0) ? -1 : 0)) * sqrt(DoF2) / dd4hep::c_light;
1057 if (level() >= algorithms::LogLevel::kDebug) {
1058 debug("--------------------");
1059 printSubHitList(input, subHitList);
1060 debug(" =");
1061
1062 Position locPos(lpos[0], lpos[1], lpos[2]);
1063 Position globPos = refVol.nominal().localToWorld(locPos);
1064 debug(" position = ({:7.2f},{:7.2f},{:7.2f}) [mm]", globPos.X() / mm, globPos.Y() / mm,
1065 globPos.Z() / mm);
1066 debug(" edep = {:.0f} [eV]", eDep / eV);
1067 debug(" time = {:.2f} [ns]", time);
1068 }
1069 return true;
1070 }
1071 void MPGDTrackerDigi::printSubHitList(const Input& input, std::vector<int>& subHitList) const {
1072 const auto [headers, sim_hits] = input;
1073 const double edmm = edm4eic::unit::mm;
1074 using edm4eic::unit::eV, edm4eic::unit::GeV;
1075 int ldx = 0;
1076 for (int kdx : subHitList) {
1077 const edm4hep::SimTrackerHit& sim_hp = sim_hits->at(kdx);
1078 CellID cIDk = sim_hp.getCellID();
1079 CellID hIDk = cIDk >> 32, vIDk = cIDk & m_volumeBits;
1080 if (ldx == 0) {
1081 debug("Hit cellID{:d} = 0x{:08x}, 0x{:08x}", ldx++, hIDk, vIDk);
1082 } else {
1083 debug(" + cellID{:d} = 0x{:08x}, 0x{:08x}", ldx++, hIDk, vIDk);
1084 }
1085 debug(" position = ({:7.2f},{:7.2f},{:7.2f}) [mm]", sim_hp.getPosition().x / edmm,
1086 sim_hp.getPosition().y / edmm, sim_hp.getPosition().z / edmm);
1087 debug(" xy_radius = {:.2f}",
1088 std::hypot(sim_hp.getPosition().x, sim_hp.getPosition().y) / edmm);
1089 debug(" momentum = ({:.2f}, {:.2f}, {:.2f}) [GeV]", sim_hp.getMomentum().x / GeV,
1090 sim_hp.getMomentum().y / GeV, sim_hp.getMomentum().z / GeV);
1091 debug(" edep = {:.0f} [eV]", sim_hp.getEDep() / eV);
1092 debug(" time = {:.2f} [ns]", sim_hp.getTime());
1093 }
1094 }
1095
1096 void getLocalPosMom(const edm4hep::SimTrackerHit& sim_hit, const TGeoHMatrix& toModule,
1097 double* lpos, double* lmom) {
1098 const edm4hep::Vector3d& pos = sim_hit.getPosition();
1099
1100 const double edmm = edm4eic::unit::mm, ed2dd = dd4hep::mm / edmm;
1101 const double gpos[3] = {pos.x * ed2dd, pos.y * ed2dd, pos.z * ed2dd};
1102 const edm4hep::Vector3f& mom = sim_hit.getMomentum();
1103 const double gmom[3] = {mom.x, mom.y, mom.z};
1104 toModule.MasterToLocal(gpos, lpos);
1105 toModule.MasterToLocalVect(gmom, lmom);
1106 }
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125 unsigned int MPGDTrackerDigi::cTraversing(const double* lpos, const double* lmom, double path,
1126 bool isSecondary,
1127 double rMin, double rMax,
1128 double dZ, double startPhi,
1129 double endPhi,
1130 double lintos[][3], double louts[][3], double* lpini,
1131 double* lpend) const {
1132 unsigned int status = 0;
1133 double Mx = lpos[0], My = lpos[1], Mz = lpos[2], M2 = Mx * Mx + My * My;
1134 double Px = lmom[0], Py = lmom[1], Pz = lmom[2];
1135
1136 double tIn = 0, tOut = 0;
1137 for (double phi : {startPhi, endPhi}) {
1138
1139 double Ux = cos(phi), Uy = sin(phi);
1140 double D = Px * Uy - Py * Ux;
1141 if (D) {
1142 double t = (My * Ux - Mx * Uy) / D;
1143 double Ex = Mx + t * Px, Ey = My + t * Py, Ez = Mz + t * Pz;
1144 double rE = sqrt(Ex * Ex + Ey * Ey), phiE = atan2(Ey, Ex);
1145
1146
1147 if (rMin < rE && rE < rMax && fabs(Ez) < dZ && fabs(phiE - phi) < 1) {
1148 if (t < 0) {
1149 status |= 0x10;
1150 tIn = t;
1151 } else {
1152 status |= 0x20;
1153 tOut = t;
1154 }
1155 }
1156 }
1157 }
1158
1159 double zLow = -dZ, zUp = +dZ;
1160 for (double Z : {zLow, zUp}) {
1161
1162 if (Pz) {
1163 double t = (Z - Mz) / Pz;
1164 double Ex = Mx + t * Px, Ey = My + t * Py, rE = sqrt(Ex * Ex + Ey * Ey);
1165 double phi = atan2(Ey, Ex);
1166 if (rMin < rE && rE < rMax && startPhi < phi && phi < endPhi) {
1167 if (t < 0) {
1168 if (!(status & 0x10) || ((status & 0x10) && t > tIn)) {
1169 status |= 0x10;
1170 tIn = t;
1171 }
1172 } else if (t > 0) {
1173 if (!(status & 0x20) || ((status & 0x20) && t < tOut)) {
1174 status |= 0x20;
1175 tOut = t;
1176 }
1177 }
1178 }
1179 }
1180 }
1181
1182 double ts[3 ][2 ] = {
1183 {0, 0}, {0, 0}, {tIn, tOut}};
1184 double a = Px * Px + Py * Py, b = Px * Mx + Py * My;
1185 for (int lu = 0; lu < 2; lu++) {
1186 double R;
1187 unsigned int statGene;
1188 if (lu == 1) {
1189 R = rMax;
1190 statGene = 0x4;
1191 } else {
1192 R = rMin;
1193 statGene = 0x1;
1194 }
1195 double c = M2 - R * R;
1196 if (!a) {
1197 if ((status & 0x30) != 0x30)
1198 status |= 0x1000;
1199 continue;
1200 } else if (!c) {
1201 status |= 0x2000;
1202 continue;
1203 } else {
1204 double det = b * b - a * c;
1205 if (det < 0) {
1206 if (lu == 1) {
1207 status |= 0x4000;
1208 continue;
1209 }
1210 } else {
1211 double sqdet = sqrt(det);
1212 for (int is = 0; is < 2; is++) {
1213 int s = 1 - 2 * is;
1214 double t = (-b + s * sqdet) / a;
1215 double Ix = Mx + t * Px, Iy = My + t * Py, Iz = Mz + t * Pz, phi = atan2(Iy, Ix);
1216 if (fabs(Iz) > dZ || phi < startPhi || endPhi < phi)
1217 continue;
1218 if (t < 0) {
1219
1220
1221
1222
1223 if (status & statGene) {
1224 double tPrv = ts[lu][0];
1225 if (t > tPrv) {
1226 ts[lu][0] = t;
1227 ts[lu][1] = tPrv;
1228 }
1229 status |= statGene << 1;
1230 } else {
1231 ts[lu][0] = t;
1232 status |= statGene;
1233 }
1234 } else {
1235 if (status & statGene << 1) {
1236 double tPrv = ts[lu][1];
1237 if (t < tPrv) {
1238 ts[lu][1] = t;
1239 ts[lu][0] = tPrv;
1240 }
1241 status |= statGene;
1242 } else {
1243 ts[lu][1] = t;
1244 status |= statGene << 1;
1245 }
1246 }
1247 }
1248 }
1249 }
1250 }
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268 bool canReEnter = (status & 0x3) == 0x3;
1269 double norm = sqrt(a + Pz * Pz);
1270 if (canReEnter) {
1271 bool doesReEnter = true;
1272 for (int i12 = 0; i12 < 2; i12++) {
1273 double t = ts[0][i12];
1274 int s = t > 0 ? +1 : -1;
1275 const double tolerance = 20 * dd4hep::um;
1276 if (path / 2 - s * t * norm < tolerance)
1277 doesReEnter = false;
1278 }
1279 if (doesReEnter) {
1280 status &= ~0x3;
1281 canReEnter = false;
1282 }
1283 }
1284 double rHit = sqrt(M2);
1285 if (rHit < rMin && !canReEnter) {
1286 unsigned int statvs = status;
1287 if (statvs & 0x1) {
1288 status &= ~0x1;
1289 status |= 0x2;
1290 ts[0][1] = -ts[0][0];
1291 }
1292 if (statvs & 0x2) {
1293 status &= ~0x2;
1294 status |= 0x1;
1295 ts[0][0] = -ts[0][1];
1296 }
1297 } else if (rHit > rMax) {
1298 unsigned int statvs = status;
1299 if (statvs & 0x4) {
1300 status &= ~0x4;
1301 status |= 0x8;
1302 ts[1][1] = -ts[1][0];
1303 }
1304 if (statvs & 0x8) {
1305 status &= ~0x8;
1306 status |= 0x4;
1307 ts[1][0] = -ts[1][1];
1308 }
1309 }
1310 for (int lu = 0; lu < 2; lu++) {
1311 unsigned int statGene = lu ? 0x4 : 0x1;
1312 if (status & statGene) {
1313 double t = ts[lu][0];
1314 if (t < 0) {
1315 if (status & 0x10) {
1316 if (lu == 1 || !canReEnter) {
1317 status |= 0x10000;
1318 status &= ~statGene;
1319 } else {
1320 if (t < tIn)
1321 status |= 0x10000;
1322 }
1323 }
1324 } else {
1325 if (status & 0x20) {
1326 if (lu == 1 || !canReEnter) {
1327 status |= 0x20000;
1328 status &= ~statGene;
1329 } else {
1330 if (t > tOut)
1331 status |= 0x20000;
1332 }
1333 }
1334 }
1335 }
1336 if (status & statGene << 1) {
1337 double t = ts[lu][1];
1338 if (t < 0) {
1339 if (status & 0x10) {
1340 if (lu == 1 || !canReEnter) {
1341 status |= 0x40000;
1342 status &= ~(statGene << 1);
1343 } else {
1344 if (t < tIn)
1345 status |= 0x40000;
1346 }
1347 }
1348 } else {
1349 if (status & 0x20) {
1350 if (lu == 1 || !canReEnter) {
1351 status |= 0x80000;
1352 status &= ~(statGene << 1);
1353 } else {
1354 if (t > tOut)
1355 status |= 0x80000;
1356 }
1357 }
1358 }
1359 }
1360 }
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370 if (canReEnter)
1371 status |= 0x100;
1372 double at = path / 2 / norm;
1373 unsigned int statws = 0;
1374 for (int is = 0; is < 2; is++) {
1375 int s = 1 - 2 * is;
1376 double Ix = s * at * Px, Iy = s * at * Py, Iz = s * at * Pz;
1377 for (int lu = 0; lu < 2; lu++) {
1378 unsigned int statvs = lu ? 0x4 : 0x1;
1379 for (int io = 0; io < 2; io++) {
1380 statvs <<= io;
1381 if (status & statvs) {
1382 double t = ts[lu][io];
1383 if (t * s < 0)
1384 continue;
1385 double dIx = t * Px - Ix, dIy = t * Py - Iy, dIz = t * Pz - Iz;
1386 double dist = sqrt(dIx * dIx + dIy * dIy + dIz * dIz);
1387
1388 double tolerance = m_toleranceFactor(norm) * 20 * dd4hep::um;
1389 if (dist < tolerance)
1390 statws |= statvs;
1391 }
1392 }
1393 }
1394 }
1395 if (!(statws & 0x5))
1396 status &= ~0x5;
1397 if (!(statws & 0xa))
1398 status &= ~0xa;
1399
1400
1401
1402 if (((status & 0x5) == 0x1 || (status & 0x5) == 0x4) && !isSecondary) {
1403 double tIn = (status & 0x1) ? ts[0][0] : ts[1][0];
1404 lpini[0] = Mx + tIn * Px;
1405 lpini[1] = My + tIn * Py;
1406 lpini[2] = Mz + tIn * Pz;
1407 } else {
1408 lpini[0] = Mx - at * Px;
1409 lpini[1] = My - at * Py;
1410 lpini[2] = Mz - at * Pz;
1411 }
1412 if (((status & 0xa) == 0x2 || (status & 0xa) == 0x8) && !isSecondary) {
1413 double tOut = (status & 0x2) ? ts[0][1] : ts[1][1];
1414 lpend[0] = Mx + tOut * Px;
1415 lpend[1] = My + tOut * Py;
1416 lpend[2] = Mz + tOut * Pz;
1417 } else {
1418 lpend[0] = Mx + at * Px;
1419 lpend[1] = My + at * Py;
1420 lpend[2] = Mz + at * Pz;
1421 }
1422
1423 for (int lu = 0; lu < 2; lu++) {
1424 unsigned int statvs = lu ? 0x4 : 0x1;
1425 double tIn = ts[lu][0], tOut = ts[lu][1];
1426 if (status & statvs) {
1427 lintos[lu][0] = Mx + tIn * Px;
1428 lintos[lu][1] = My + tIn * Py;
1429 lintos[lu][2] = Mz + tIn * Pz;
1430 }
1431 statvs <<= 1;
1432 if (status & statvs) {
1433 louts[lu][0] = Mx + tOut * Px;
1434 louts[lu][1] = My + tOut * Py;
1435 louts[lu][2] = Mz + tOut * Pz;
1436 }
1437 }
1438 return status;
1439 }
1440 unsigned int MPGDTrackerDigi::bTraversing(const double* lpos, const double* lmom, double ref2Cur,
1441 double path,
1442 bool isSecondary,
1443 double dZ,
1444 double dX, double dY,
1445 double lintos[][3], double louts[][3], double* lpini,
1446 double* lpend) const {
1447 unsigned int status = 0;
1448 double Mx = lpos[0], My = lpos[1], Mxy[2] = {Mx, My};
1449 double Px = lmom[0], Py = lmom[1], Pxy[2] = {Px, Py};
1450 double Mz = lpos[2] + ref2Cur, Pz = lmom[2];
1451
1452 double tIn = 0, tOut = 0;
1453 double xyLow[2] = {-dX, -dY}, xyUp[2] = {+dX, +dY};
1454 for (int xy = 0; xy < 2; xy++) {
1455 int yx = 1 - xy;
1456 double a_Low = xyLow[xy], a_Up = xyUp[xy], Ma = Mxy[xy], Pa = Pxy[xy];
1457 double b_Low = xyLow[yx], b_Up = xyUp[yx], Mb = Mxy[yx], Pb = Pxy[yx];
1458 for (double A : {a_Low, a_Up}) {
1459
1460 if (Pa) {
1461 double t = (A - Ma) / Pa;
1462 double Eb = Mb + t * Pb, Ez = Mz + t * Pz;
1463 if (b_Low < Eb && Eb < b_Up && fabs(Ez) < dZ) {
1464 if (t < 0) {
1465 if (!(status & 0x10) || ((status & 0x10) && t > tIn)) {
1466 status |= 0x10;
1467 tIn = t;
1468 }
1469 } else if (t > 0) {
1470 if (!(status & 0x20) || ((status & 0x20) && t < tOut)) {
1471 status |= 0x20;
1472 tOut = t;
1473 }
1474 }
1475 }
1476 }
1477 }
1478 }
1479
1480 for (int lu = 0; lu < 2; lu++) {
1481 int s = 2 * lu - 1;
1482 double Z = s * dZ;
1483 unsigned int statGene = lu ? 0x4 : 0x1;
1484
1485 if (Pz) {
1486 double t = (Z - Mz) / Pz;
1487 if (t < 0) {
1488 if (!(status & 0x10) || ((status & 0x10) && t > tIn)) {
1489 status |= statGene;
1490 tIn = t;
1491 }
1492 } else if (t > 0) {
1493 if (!(status & 0x20) || ((status & 0x20) && t < tOut)) {
1494 status |= statGene << 1;
1495 tOut = t;
1496 }
1497 }
1498 }
1499 }
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509 double norm = sqrt(Px * Px + Py * Py + Pz * Pz), at = path / 2 / norm;
1510 unsigned int statws = 0;
1511 for (int is = 0; is < 2; is++) {
1512 int s = 1 - 2 * is;
1513 double Ix = s * at * Px, Iy = s * at * Py, Iz = s * at * Pz;
1514 for (int lu = 0; lu < 2; lu++) {
1515 unsigned int statvs = lu ? 0x4 : 0x1;
1516 for (int io = 0; io < 2; io++) {
1517 statvs <<= io;
1518 if (status & statvs) {
1519 double t = io ? tOut : tIn;
1520 if (t * s < 0)
1521 continue;
1522 double dIx = t * Px - Ix, dIy = t * Py - Iy, dIz = t * Pz - Iz;
1523 double dist = sqrt(dIx * dIx + dIy * dIy + dIz * dIz);
1524
1525 double tolerance = m_toleranceFactor(norm) * 20 * dd4hep::um;
1526 if (dist < tolerance)
1527 statws |= statvs;
1528 }
1529 }
1530 }
1531 }
1532 if (!(statws & 0x5))
1533 status &= ~0x5;
1534 if (!(statws & 0xa))
1535 status &= ~0xa;
1536
1537 Mz -= ref2Cur;
1538
1539
1540 if ((status & 0x5) && !isSecondary) {
1541 lpini[0] = Mx + tIn * Px;
1542 lpini[1] = My + tIn * Py;
1543 lpini[2] = Mz + tIn * Pz;
1544 } else {
1545 lpini[0] = Mx - at * Px;
1546 lpini[1] = My - at * Py;
1547 lpini[2] = Mz - at * Pz;
1548 }
1549 if ((status & 0xa) && !isSecondary) {
1550 lpend[0] = Mx + tOut * Px;
1551 lpend[1] = My + tOut * Py;
1552 lpend[2] = Mz + tOut * Pz;
1553 } else {
1554 lpend[0] = Mx + at * Px;
1555 lpend[1] = My + at * Py;
1556 lpend[2] = Mz + at * Pz;
1557 }
1558
1559 for (int lu = 0; lu < 2; lu++) {
1560 unsigned int statvs = lu ? 0x4 : 0x1;
1561 if (status & statvs) {
1562 lintos[lu][0] = Mx + tIn * Px;
1563 lintos[lu][1] = My + tIn * Py;
1564 lintos[lu][2] = Mz + tIn * Pz;
1565 }
1566 statvs <<= 1;
1567 if (status & statvs) {
1568 louts[lu][0] = Mx + tOut * Px;
1569 louts[lu][1] = My + tOut * Py;
1570 louts[lu][2] = Mz + tOut * Pz;
1571 }
1572 }
1573 return status;
1574 }
1575
1576
1577 bool cExtrapolate(const double* lpos, const double* lmom,
1578 double rT,
1579 double* lext)
1580 {
1581 bool ok = false;
1582 double Mx = lpos[0], My = lpos[1], Mz = lpos[2], M2 = Mx * Mx + My * My;
1583 double Px = lmom[0], Py = lmom[1], Pz = lmom[2];
1584 double a = Px * Px + Py * Py, b = Px * Mx + Py * My, c = M2 - rT * rT;
1585 double tF = 0;
1586 if (!c)
1587 ok = true;
1588 else if (a) {
1589 double det = b * b - a * c;
1590 if (det >= 0) {
1591 double sqdet = sqrt(det);
1592 for (int is = 0; is < 2; is++) {
1593 int s = 1 - 2 * is;
1594 double t = (-b + s * sqdet) / a, norm = sqrt(a + Pz * Pz);
1595
1596 if (t * norm < -dd4hep::nm)
1597 continue;
1598 if (!ok ||
1599
1600 (ok && fabs(t) < fabs(tF))) {
1601 tF = t;
1602 ok = true;
1603 }
1604 }
1605 }
1606 }
1607 if (ok) {
1608 lext[0] = Mx + tF * Px;
1609 lext[1] = My + tF * Py;
1610 lext[2] = Mz + tF * Pz;
1611 }
1612 return ok;
1613 }
1614 bool bExtrapolate(const double* lpos, const double* lmom,
1615 double zT,
1616 double* lext)
1617 {
1618 bool ok = false;
1619 double Mx = lpos[0], My = lpos[1], Mz = lpos[2];
1620 double Px = lmom[0], Py = lmom[1], Pz = lmom[2], norm = sqrt(Px * Px + Py * Py + Pz * Pz);
1621 double tF = 0;
1622 if (Pz) {
1623 tF = (zT - Mz) / Pz;
1624
1625 ok = tF * norm > -dd4hep::nm;
1626 }
1627 if (ok) {
1628 lext[0] = Mx + tF * Px;
1629 lext[1] = My + tF * Py;
1630 lext[2] = Mz + tF * Pz;
1631 }
1632 return ok;
1633 }
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644 unsigned int MPGDTrackerDigi::cExtension(double const* lpos, double const* lmom,
1645 double rT, int direction,
1646 double dZ, double startPhi,
1647 double endPhi,
1648 double* lext) const {
1649 unsigned int status = 0;
1650 double Mx = lpos[0], My = lpos[1], Mz = lpos[2];
1651 double Px = lmom[0], Py = lmom[1], Pz = lmom[2], norm = sqrt(Px * Px + Py * Py + Pz * Pz);
1652
1653
1654 const double margin = 10 * dd4hep::um;
1655 double t = direction * margin / norm;
1656 Mx += t * Px;
1657 My += t * Py;
1658 Mz += t * Pz;
1659 double M2 = Mx * Mx + My * My, rIni = sqrt(M2), rLow, rUp;
1660 if (rIni < rT) {
1661 rLow = rIni;
1662 rUp = rT;
1663 } else {
1664 rLow = rT;
1665 rUp = rIni;
1666 }
1667
1668 double tF = 0;
1669 for (double phi : {startPhi, endPhi}) {
1670
1671 double Ux = cos(phi), Uy = sin(phi);
1672 double D = Px * Uy - Py * Ux;
1673 if (D) {
1674 double t = (My * Ux - Mx * Uy) / D;
1675 if (t * direction < 0)
1676 continue;
1677 double Ex = Mx + t * Px, Ey = My + t * Py, Ez = Mz + t * Pz;
1678 double rE = sqrt(Ex * Ex + Ey * Ey), phiE = atan2(Ey, Ex);
1679
1680 if (rLow < rE && rE < rUp && fabs(Ez) < dZ && fabs(phiE - phi) < 1) {
1681 status |= 0x1;
1682 tF = t;
1683 }
1684 }
1685 }
1686
1687 double zLow = -dZ, zUp = +dZ;
1688 for (double Z : {zLow, zUp}) {
1689
1690 if (Pz) {
1691 double t = (Z - Mz) / Pz;
1692 if (t * direction < 0)
1693 continue;
1694 double Ex = Mx + t * Px, Ey = My + t * Py, rE = sqrt(Ex * Ex + Ey * Ey);
1695 double phi = atan2(Ey, Ex);
1696 if (rLow < rE && rE < rUp && startPhi < phi && phi < endPhi) {
1697 if (t < 0) {
1698 if (!status || (status && t > tF)) {
1699 status |= 0x1;
1700 tF = t;
1701 }
1702 } else if (t > 0) {
1703 if (!status || (status && t < tF)) {
1704 status |= 0x1;
1705 tF = t;
1706 }
1707 }
1708 }
1709 }
1710 }
1711
1712 if (!status) {
1713 double a = Px * Px + Py * Py, b = Px * Mx + Py * My, c = M2 - rT * rT;
1714 if (!a) {
1715 status |= 0x1000;
1716 } else if (!c) {
1717 status |= 0x2000;
1718 } else {
1719 double det = b * b - a * c;
1720 if (det >= 0) {
1721 double sqdet = sqrt(det);
1722 for (int is = 0; is < 2; is++) {
1723 int s = 1 - 2 * is;
1724 double t = (-b + s * sqdet) / a;
1725 if (t * direction < 0)
1726 continue;
1727 double Ix = Mx + t * Px, Iy = My + t * Py, Iz = Mz + t * Pz, phi = atan2(Iy, Ix);
1728 if (fabs(Iz) > dZ || phi < startPhi || endPhi < phi)
1729 continue;
1730 if (!(status & 0x1) ||
1731
1732 ((status & 0x1) && fabs(t) < fabs(tF))) {
1733 tF = t;
1734 status |= 0x1;
1735 }
1736 }
1737 }
1738 }
1739 }
1740 if (status & 0x1) {
1741 lext[0] = Mx + tF * Px;
1742 lext[1] = My + tF * Py;
1743 lext[2] = Mz + tF * Pz;
1744 }
1745 return status;
1746 }
1747 unsigned int MPGDTrackerDigi::bExtension(const double* lpos, const double* lmom,
1748 double zT, int direction,
1749 double dX, double dY,
1750 double* lext) const {
1751 unsigned int status = 0;
1752 double Mx = lpos[0], My = lpos[1], Mxy[2] = {Mx, My};
1753 double Px = lmom[0], Py = lmom[1], Pxy[2] = {Px, Py};
1754 double Mz = lpos[2], Pz = lmom[2];
1755 double norm = sqrt(Px * Px + Py * Py + Pz * Pz);
1756
1757
1758 const double margin = 10 * dd4hep::um;
1759 double t = direction * margin / norm;
1760 Mx += t * Px;
1761 My += t * Py;
1762 Mz += t * Pz;
1763 double &zIni = Mz, zLow, zUp;
1764 if (zIni < zT) {
1765 zLow = zIni;
1766 zUp = zT;
1767 } else {
1768 zLow = zT;
1769 zUp = zIni;
1770 }
1771
1772 double tF = 0;
1773 double xyLow[2] = {-dX, -dY}, xyUp[2] = {+dX, +dY};
1774 for (int xy = 0; xy < 2; xy++) {
1775 int yx = 1 - xy;
1776 double a_Low = xyLow[xy], a_Up = xyUp[xy], Ma = Mxy[xy], Pa = Pxy[xy];
1777 double b_Low = xyLow[yx], b_Up = xyUp[yx], Mb = Mxy[yx], Pb = Pxy[yx];
1778 for (double A : {a_Low, a_Up}) {
1779
1780 if (Pa) {
1781 double t = (A - Ma) / Pa;
1782 if (t * direction < 0)
1783 continue;
1784 double Eb = Mb + t * Pb, Ez = Mz + t * Pz;
1785 if (zLow < Ez && Ez < zUp && b_Low < Eb && Eb < b_Up) {
1786 if (!status || (status && fabs(t) < fabs(tF))) {
1787 status |= 0x1;
1788 tF = t;
1789 }
1790 }
1791 }
1792 }
1793 }
1794
1795 if (!status) {
1796 if (Pz) {
1797 tF = (zT - Mz) / Pz;
1798 if (tF * direction > 0)
1799 status = 0x1;
1800 }
1801 }
1802 if (status) {
1803 lext[0] = Mx + tF * Px;
1804 lext[1] = My + tF * Py;
1805 lext[2] = Mz + tF * Pz;
1806 }
1807 return status;
1808 }
1809
1810 double getRef2Cur(DetElement refVol, DetElement curVol) {
1811
1812
1813 const TGeoHMatrix toRefVol = refVol.nominal().worldTransformation();
1814 const TGeoHMatrix toCurVol = curVol.nominal().worldTransformation();
1815 const double* TRef = toRefVol.GetTranslation();
1816 const double* TCur = toCurVol.GetTranslation();
1817
1818 double gdT[3];
1819 for (int i = 0; i < 3; i++)
1820 gdT[i] = TRef[i] - TCur[i];
1821 double ldT[3];
1822 toRefVol.MasterToLocalVect(gdT, ldT);
1823 return ldT[2];
1824 }
1825
1826 std::string inconsistency(const edm4hep::EventHeader& event, unsigned int status, CellID cID,
1827 const double* lpos, const double* lmom) {
1828 using edm4eic::unit::GeV, dd4hep::mm;
1829 return fmt::format("Event {}#{}, SimHit 0x{:016x} @ {:.2f},{:.2f},{:.2f} mm, P = "
1830 "{:.2f},{:.2f},{:.2f} GeV inconsistency 0x{:x}",
1831 event.getRunNumber(), event.getEventNumber(), cID, lpos[0] / mm, lpos[1] / mm,
1832 lpos[2] / mm, lmom[0] / GeV, lmom[1] / GeV, lmom[2] / GeV, status);
1833 }
1834 std::string oddity(const edm4hep::EventHeader& event, unsigned int status, double dist, CellID cID,
1835 const double* lpos, const double* lmom, CellID cJD, const double* lpoj,
1836 const double* lmoj) {
1837 using edm4eic::unit::GeV, dd4hep::mm;
1838 return fmt::format("Event {}#{}, Bizarre SimHit sequence: 0x{:016x} @ {:.4f},{:.4f},{:.4f} mm, P "
1839 "= {:.2f},{:.2f},{:.2f} GeV and 0x{:016x} @ {:.4f},{:.4f},{:.4f} mm, P = "
1840 "{:.2f},{:.2f},{:.2f} GeV: status 0x{:x}, distance {:.4f}",
1841 event.getRunNumber(), event.getEventNumber(), cID, lpos[0] / mm, lpos[1] / mm,
1842 lpos[2] / mm, lmom[0] / GeV, lmom[1] / GeV, lmom[2] / GeV, cJD, lpoj[0] / mm,
1843 lpoj[1] / mm, lpoj[2] / mm, lmoj[0] / GeV, lmoj[1] / GeV, lmoj[2] / GeV,
1844 status, dist);
1845 }
1846
1847 bool MPGDTrackerDigi::samePMO(const edm4hep::SimTrackerHit& sim_hit,
1848 const edm4hep::SimTrackerHit& sim_hjt) const {
1849
1850
1851
1852
1853 #if EDM4HEP_BUILD_VERSION >= EDM4HEP_VERSION(0, 99, 0)
1854 bool sameParticle = sim_hjt.getParticle() == sim_hit.getParticle();
1855 #else
1856 bool sameParticle = sim_hjt.getMCParticle() == sim_hit.getMCParticle();
1857 #endif
1858
1859 CellID vID = sim_hit.getCellID() & m_volumeBits;
1860 CellID refID = vID & m_moduleBits;
1861 CellID vJD = sim_hjt.getCellID() & m_volumeBits;
1862 CellID refJD = vJD & m_moduleBits;
1863 bool sameModule = refJD == refID;
1864
1865
1866
1867 bool isSecondary = sim_hit.isProducedBySecondary();
1868 bool jsSecondary = sim_hjt.isProducedBySecondary();
1869 bool sameOrigin = jsSecondary == isSecondary;
1870 return sameParticle && sameModule && sameOrigin;
1871 }
1872
1873 double outInDistance(int shape, int orientation, double lintos[][3], double louts[][3],
1874 double* lmom, double* lmoj) {
1875
1876 bool ok;
1877 double lExt[3];
1878 double lmOI[3];
1879 for (int i = 0; i < 3; i++)
1880 lmOI[i] = (lmom[i] + lmoj[i]) / 2;
1881 double *lOut, *lInto;
1882 if (orientation > 0) {
1883 lOut = louts[1];
1884 lInto = lintos[0];
1885 } else if (orientation < 0) {
1886 lOut = louts[0];
1887 lInto = lintos[1];
1888 } else {
1889 lOut = louts[0];
1890 lInto = lintos[0];
1891 }
1892 if (shape == 0) {
1893 double rInto = sqrt(lInto[0] * lInto[0] + lInto[1] * lInto[1]);
1894 ok = cExtrapolate(lOut, lmOI, rInto, lExt);
1895 } else {
1896 ok = bExtrapolate(lOut, lmOI, lInto[2], lExt);
1897 }
1898 if (ok) {
1899 double dist2 = 0;
1900 for (int i = 0; i < 3; i++) {
1901 double d = lExt[i] - lInto[i];
1902 dist2 += d * d;
1903 }
1904 return sqrt(dist2);
1905 } else
1906 return -1;
1907 }
1908
1909 unsigned int MPGDTrackerDigi::extendHit(CellID refID, std::vector<std::uint64_t>& cIDs,
1910 int direction, double* lpini, double* lmini, double* lpend,
1911 double* lmend) const {
1912 unsigned int status = 0;
1913 const VolumeManager& volman = m_detector->volumeManager();
1914 DetElement refVol = volman.lookupDetElement(refID);
1915 const auto& shape = refVol.solid();
1916 double *lpoE, *lmoE;
1917 if (direction < 0) {
1918 lpoE = lpini;
1919 lmoE = lmini;
1920 } else {
1921 lpoE = lpend;
1922 lmoE = lmend;
1923 }
1924
1925
1926
1927
1928
1929
1930
1931 for (int rankE : {0, 4}) {
1932 CellID vIDE = refID | m_stripIDs[rankE];
1933 int alreadyThere = 0;
1934 for (int i = 0; i < (int)cIDs.size(); i++) {
1935 if ((cIDs[i] & m_volumeBits) == vIDE) {
1936 alreadyThere = 1;
1937 break;
1938 }
1939 }
1940 if (alreadyThere)
1941 continue;
1942 if (std::find(cIDs.begin(), cIDs.end(), vIDE) != cIDs.end())
1943 continue;
1944 DetElement volE = volman.lookupDetElement(vIDE);
1945 double lext[3];
1946 if (std::string_view{shape.type()} == "TGeoTubeSeg") {
1947 const Tube& tExt = volE.solid();
1948 double R = rankE == 0 ? tExt.rMin() : tExt.rMax();
1949 double startPhi = tExt.startPhi() * radian;
1950 startPhi -= 2 * TMath::Pi();
1951 double endPhi = tExt.endPhi() * radian;
1952 endPhi -= 2 * TMath::Pi();
1953 double dZ = tExt.dZ();
1954 status = cExtension(lpoE, lmoE, R, direction, dZ, startPhi, endPhi, lext);
1955 } else if (std::string_view{shape.type()} == "TGeoBBox") {
1956 double ref2E = getRef2Cur(refVol, volE);
1957 const Box& bExt = volE.solid();
1958 double Z = rankE == 0 ? -bExt.z() : +bExt.z();
1959 Z -= ref2E;
1960 double dX = bExt.x(), dY = bExt.y();
1961 status = bExtension(lpoE, lmoE, Z, direction, dX, dY, lext);
1962 } else {
1963 critical(R"(Bad input data: CellID {:x} has invalid shape "{}")", refID, shape.type());
1964 throw std::runtime_error(R"(Inconsistency: Inappropriate SimHits fed to "MPGDTrackerDigi".)");
1965 }
1966 if (status != 0x1)
1967 continue;
1968 if (direction < 0) {
1969 for (int i = 0; i < 3; i++)
1970 lpini[i] = lext[i];
1971 } else {
1972 for (int i = 0; i < 3; i++)
1973 lpend[i] = lext[i];
1974 }
1975 break;
1976 }
1977 return status;
1978 }
1979
1980 bool MPGDTrackerDigi::denyExtension(const edm4hep::SimTrackerHit& sim_hit, double depth) const {
1981
1982
1983
1984 CellID vID = sim_hit.getCellID() & m_volumeBits;
1985 bool isHelperVolume = m_stripRank(vID) != 0 && m_stripRank(vID) != 4;
1986
1987
1988 const double fraction = .10;
1989 const double edmm = edm4eic::unit::mm, ed2dd = dd4hep::mm / edmm;
1990 bool smallPathLength = sim_hit.getPathLength() * ed2dd < fraction * depth;
1991 return isHelperVolume || smallPathLength;
1992 }
1993
1994 void MPGDTrackerDigi::flagUnexpected(const edm4hep::EventHeader& event, int shape, double expected,
1995 const edm4hep::SimTrackerHit& sim_hit, double* lpini,
1996 double* lpend, double* lpos, double* lmom) const {
1997
1998
1999
2000
2001
2002
2003
2004
2005
2006
2007 double Rnew2 = 0, Znew, diff2 = 0;
2008 for (int i = 0; i < 3; i++) {
2009 double neu = (lpini[i] + lpend[i]) / 2, alt = lpos[i];
2010 double d = neu - alt;
2011 diff2 += d * d;
2012 if (i != 2)
2013 Rnew2 += neu * neu;
2014 if (i == 2)
2015 Znew = neu;
2016 }
2017 double found = shape ? Znew : sqrt(Rnew2), residual = found - expected;
2018 bool isSecondary = sim_hit.isProducedBySecondary();
2019 bool isPrimary =
2020 !isSecondary && sqrt(lmom[0] * lmom[0] + lmom[1] * lmom[1] + lmom[2] * lmom[2]) > .1 * GeV;
2021 if ((fabs(residual) > .000001 && isPrimary) || (sqrt(diff2) > .000001 && isSecondary)) {
2022 debug("Event {}#{}, SimHit 0x{:016x} origin {:d}: d{:c} = {:.5f} diff = {:.5f}",
2023 event.getRunNumber(), event.getEventNumber(), sim_hit.getCellID(), isSecondary,
2024 shape ? 'Z' : 'R', residual, sqrt(diff2));
2025 }
2026 }
2027
2028
2029
2030
2031
2032 int MPGDTrackerDigi::get2HitCluster(CellID refID,
2033 Position& locPos,
2034 double surfPos[2],
2035 int pn,
2036 std::default_random_engine& generator, Cluster& cluster) const {
2037
2038
2039 CellID stripID = m_stripIDs[pn ? 3 : 1];
2040 const Position dummy(0, 0, 0);
2041 CellID masterID = m_seg->cellID(locPos, dummy, refID | stripID);
2042
2043 const StripParameters* pars;
2044 CellID sensorStripID = masterID & m_sensorStripBits;
2045 try {
2046 pars = &m_stripParameters.at(sensorStripID);
2047 } catch (const std::out_of_range& oor) {
2048 critical(R"(Error retrieving StripParameters for readout "{}", cellID 0x{:0>16x}: {}.)",
2049 m_cfg.readout, masterID, oor.what());
2050 throw std::runtime_error("Error retrieving StripParameters");
2051 }
2052
2053
2054
2055
2056
2057 const double& sigma = m_cfg.stripResolutions[pn];
2058 const double min = pars->min, max = pars->max, pitch = pars->pitch;
2059 double hA = surfPos[pn];
2060 if (hA < min - (m_truncation - .1 * dd4hep::cm) * sigma ||
2061 hA > max + (m_truncation - .1 * dd4hep::cm) * sigma) {
2062
2063
2064
2065
2066
2067
2068 return 1;
2069 }
2070 FieldID stripNum = m_id_dec->get(masterID, pars->index);
2071 double mA = m_binToPosition(stripNum, pitch, pars->offset);
2072 std::normal_distribution<double> gaussian;
2073 auto stripGauss = [&](double abscissa, double sigma, double min, double max) {
2074
2075 double low = abscissa - m_truncation * sigma;
2076 double up = abscissa + m_truncation * sigma;
2077 if (min > low)
2078 low = min;
2079 if (max < up)
2080 up = max;
2081 double x;
2082 do {
2083 x = abscissa + gaussian(generator) * sigma;
2084 } while (x < low || up < x);
2085 return x;
2086 };
2087 double sA = stripGauss(hA, sigma, min, max);
2088
2089
2090 if (sA < min + pitch / 2 || sA > max - pitch / 2) {
2091
2092 cluster.push_back({masterID, 1});
2093 } else {
2094
2095 CellID neighID, inc = m_stripIncs[pn];
2096 double nA;
2097
2098
2099
2100
2101
2102
2103
2104
2105
2106 auto incrementID = [&](int pn) {
2107 CellID spectatorBits = pn == 0 ? ((CellID)0xffff) << 48 : ((CellID)0xffff) << 32;
2108 CellID iniID = masterID & spectatorBits;
2109 neighID = masterID + inc;
2110 neighID &= ~spectatorBits;
2111 neighID |= iniID;
2112 };
2113 auto decrementID = [&](int pn) {
2114 CellID spectatorBits = pn == 0 ? ((CellID)0xffff) << 48 : ((CellID)0xffff) << 32;
2115 CellID iniID = masterID & spectatorBits;
2116 neighID = masterID - inc;
2117 neighID &= ~spectatorBits;
2118 neighID |= iniID;
2119 };
2120 if (sA < mA) {
2121 decrementID(pn);
2122 nA = mA - pitch;
2123 } else {
2124 incrementID(pn);
2125 nA = mA + pitch;
2126 }
2127
2128 double fn = (sA - mA) / (nA - mA);
2129 cluster.push_back({masterID, 1 - fn});
2130 cluster.push_back({neighID, fn});
2131 }
2132 return 0;
2133 }
2134
2135 }