File indexing completed on 2023-11-02 03:11:30
0001
0002
0003
0004
0005
0006 #include "DataFormats/Math/interface/angle_units.h"
0007 #include "DataFormats/Math/interface/FastMath.h"
0008 #include "DataFormats/ForwardDetId/interface/HGCSiliconDetId.h"
0009 #include "SimG4CMS/Calo/interface/HGCalSD.h"
0010 #include "SimG4Core/Notification/interface/TrackInformation.h"
0011 #include "FWCore/Utilities/interface/Exception.h"
0012 #include "Geometry/HGCalCommonData/interface/HGCalDDDConstants.h"
0013 #include "Geometry/HGCalCommonData/interface/HGCalGeometryMode.h"
0014 #include "G4LogicalVolumeStore.hh"
0015 #include "G4LogicalVolume.hh"
0016 #include "G4Step.hh"
0017 #include "G4Track.hh"
0018 #include "G4ParticleTable.hh"
0019 #include "G4VProcess.hh"
0020 #include "G4Trap.hh"
0021
0022 #include <cmath>
0023 #include <fstream>
0024 #include <iomanip>
0025 #include <iostream>
0026 #include <memory>
0027 #include <sstream>
0028
0029
0030
0031 using namespace angle_units::operators;
0032
0033 HGCalSD::HGCalSD(const std::string& name,
0034 const HGCalDDDConstants* hgc,
0035 const SensitiveDetectorCatalog& clg,
0036 edm::ParameterSet const& p,
0037 const SimTrackManager* manager)
0038 : CaloSD(name,
0039 clg,
0040 p,
0041 manager,
0042 static_cast<float>(p.getParameter<edm::ParameterSet>("HGCSD").getParameter<double>("TimeSliceUnit")),
0043 p.getParameter<edm::ParameterSet>("HGCSD").getParameter<bool>("IgnoreTrackID"),
0044 ("Calibration" + name)),
0045 myName_(name),
0046 hgcons_(hgc),
0047 ps_(p),
0048 slopeMin_(0),
0049 levelT1_(99),
0050 levelT2_(99),
0051 useSimWt_(0),
0052 calibCells_(false),
0053 tan30deg_(std::tan(30.0 * CLHEP::deg)),
0054 cos30deg_(std::cos(30.0 * CLHEP::deg)) {
0055 numberingScheme_.reset(nullptr);
0056 guardRing_.reset(nullptr);
0057 guardRingPartial_.reset(nullptr);
0058 mouseBite_.reset(nullptr);
0059 cellOffset_.reset(nullptr);
0060
0061 edm::ParameterSet m_HGC = p.getParameter<edm::ParameterSet>("HGCSD");
0062 eminHit_ = m_HGC.getParameter<double>("EminHit") * CLHEP::MeV;
0063 fiducialCut_ = m_HGC.getParameter<bool>("FiducialCut");
0064 storeAllG4Hits_ = m_HGC.getParameter<bool>("StoreAllG4Hits");
0065 rejectMB_ = m_HGC.getParameter<bool>("RejectMouseBite");
0066 waferRot_ = m_HGC.getParameter<bool>("RotatedWafer");
0067 cornerMinMask_ = m_HGC.getParameter<int>("CornerMinMask");
0068 nHC_ = m_HGC.getParameter<int>("HitCollection");
0069 angles_ = m_HGC.getUntrackedParameter<std::vector<double>>("WaferAngles");
0070 missingFile_ = m_HGC.getUntrackedParameter<std::string>("MissingWaferFile");
0071 checkID_ = m_HGC.getUntrackedParameter<bool>("CheckID");
0072 verbose_ = m_HGC.getUntrackedParameter<int>("Verbosity");
0073 dd4hep_ = p.getParameter<bool>("g4GeometryDD4hepSource");
0074
0075 if (storeAllG4Hits_) {
0076 setUseMap(false);
0077 setNumberCheckedHits(0);
0078 }
0079
0080
0081 mydet_ = DetId::Forward;
0082 nameX_ = "HGCal";
0083 if (myName_.find("HitsEE") != std::string::npos) {
0084 mydet_ = DetId::HGCalEE;
0085 nameX_ = "HGCalEESensitive";
0086 } else if (myName_.find("HitsHEfront") != std::string::npos) {
0087 mydet_ = DetId::HGCalHSi;
0088 nameX_ = "HGCalHESiliconSensitive";
0089 }
0090
0091 #ifdef EDM_ML_DEBUG
0092 edm::LogVerbatim("HGCSim") << "**************************************************"
0093 << "\n"
0094 << "* *"
0095 << "\n"
0096 << "* Constructing a HGCalSD with name " << name << "\n"
0097 << "* *"
0098 << "\n"
0099 << "**************************************************";
0100 #endif
0101 edm::LogVerbatim("HGCSim") << "HGCalSD:: Threshold for storing hits: " << eminHit_ << " for " << nameX_
0102 << " detector " << mydet_ << " with " << nHC_ << " hit collections";
0103 edm::LogVerbatim("HGCSim") << "Flag for storing individual Geant4 Hits " << storeAllG4Hits_;
0104 edm::LogVerbatim("HGCSim") << "Fiducial volume cut with cut from eta/phi "
0105 << "boundary " << fiducialCut_ << " at " << distanceFromEdge_;
0106 edm::LogVerbatim("HGCSim") << "Reject MosueBite Flag: " << rejectMB_ << " cuts along " << angles_.size()
0107 << " axes: " << angles_[0] << ", " << angles_[1];
0108 }
0109
0110 double HGCalSD::getEnergyDeposit(const G4Step* aStep) {
0111 double r = aStep->GetPreStepPoint()->GetPosition().perp();
0112 double z = std::abs(aStep->GetPreStepPoint()->GetPosition().z());
0113 #ifdef EDM_ML_DEBUG
0114 G4int parCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
0115 G4String parName = aStep->GetTrack()->GetDefinition()->GetParticleName();
0116 G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
0117 edm::LogVerbatim("HGCSim") << "HGCalSD: Hit from standard path from " << lv->GetName() << " for Track "
0118 << aStep->GetTrack()->GetTrackID() << " (" << parCode << ":" << parName << ") R = " << r
0119 << " Z = " << z << " slope = " << r / z << ":" << slopeMin_;
0120 #endif
0121
0122 if (r < z * slopeMin_) {
0123 #ifdef EDM_ML_DEBUG
0124 edm::LogVerbatim("HGCSim") << "HGCalSD: Fiducial Volume cut";
0125 #endif
0126 return 0.0;
0127 }
0128
0129 double wt1 = getResponseWt(aStep->GetTrack());
0130 double wt2 = aStep->GetTrack()->GetWeight();
0131 double destep = weight_ * wt1 * (aStep->GetTotalEnergyDeposit());
0132 if (wt2 > 0)
0133 destep *= wt2;
0134 #ifdef EDM_ML_DEBUG
0135 edm::LogVerbatim("HGCSim") << "HGCalSD: weights= " << weight_ << ":" << wt1 << ":" << wt2 << " Total weight "
0136 << weight_ * wt1 * wt2 << " deStep: " << aStep->GetTotalEnergyDeposit() << ":" << destep;
0137 #endif
0138 return destep;
0139 }
0140
0141 uint32_t HGCalSD::setDetUnitId(const G4Step* aStep) {
0142 const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
0143 const G4VTouchable* touch = preStepPoint->GetTouchable();
0144 fraction_ = 1.0;
0145 calibCell_ = false;
0146 int dn = touch->GetHistoryDepth();
0147
0148 #ifdef EDM_ML_DEBUG
0149 edm::LogVerbatim("HGCSim") << "DepthsTop: " << dn << ":" << levelT1_ << ":" << levelT2_;
0150 printDetectorLevels(touch);
0151 #endif
0152
0153 G4ThreeVector hitPoint = preStepPoint->GetPosition();
0154 float globalZ = touch->GetTranslation(0).z();
0155 int iz(globalZ > 0 ? 1 : -1);
0156
0157 int layer(0), moduleLev(-1), cell(-1);
0158 if (useSimWt_ > 0) {
0159 layer = touch->GetReplicaNumber(2);
0160 moduleLev = 1;
0161 } else if (touch->GetHistoryDepth() > levelT2_) {
0162 layer = touch->GetReplicaNumber(4);
0163 cell = touch->GetReplicaNumber(1);
0164 moduleLev = 3;
0165 } else {
0166 layer = touch->GetReplicaNumber(3);
0167 moduleLev = 2;
0168 }
0169 int module = touch->GetReplicaNumber(moduleLev);
0170 if (verbose_ && (cell == -1))
0171 edm::LogVerbatim("HGCSim") << "Top " << touch->GetVolume(0)->GetName() << " Module "
0172 << touch->GetVolume(moduleLev)->GetName();
0173 #ifdef EDM_ML_DEBUG
0174 edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_ << ":"
0175 << useSimWt_ << " name " << touch->GetVolume(0)->GetName() << " layer:module:cell "
0176 << layer << ":" << moduleLev << ":" << module << ":" << cell;
0177 printDetectorLevels(touch);
0178 G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
0179 edm::LogVerbatim("HGCSim") << "Depths: " << touch->GetHistoryDepth() << " name " << touch->GetVolume(0)->GetName()
0180 << ":" << touch->GetReplicaNumber(0) << " " << touch->GetVolume(1)->GetName() << ":"
0181 << touch->GetReplicaNumber(1) << " " << touch->GetVolume(2)->GetName() << ":"
0182 << touch->GetReplicaNumber(2) << " " << touch->GetVolume(3)->GetName() << ":"
0183 << touch->GetReplicaNumber(3) << " " << touch->GetVolume(4)->GetName() << ":"
0184 << touch->GetReplicaNumber(4) << " "
0185 << " layer:module:cell " << layer << ":" << module << ":" << cell << " Material "
0186 << mat->GetName() << ":" << mat->GetRadlen();
0187 #endif
0188
0189 if (aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() > 100000.)
0190 return 0;
0191
0192 uint32_t id = setDetUnitId(layer, module, cell, iz, hitPoint);
0193 if ((rejectMB_ || fiducialCut_) && id != 0) {
0194 auto uv = HGCSiliconDetId(id).waferUV();
0195 #ifdef EDM_ML_DEBUG
0196 edm::LogVerbatim("HGCSim") << "ID " << std::hex << id << std::dec << " " << HGCSiliconDetId(id);
0197 #endif
0198 G4ThreeVector local = (touch->GetHistory()->GetTransform(dn - moduleLev).TransformPoint(hitPoint));
0199 #ifdef EDM_ML_DEBUG
0200 edm::LogVerbatim("HGCSim") << "Global Point " << hitPoint << " Down0 "
0201 << touch->GetHistory()->GetTransform(dn).TransformPoint(hitPoint) << " Down1 "
0202 << touch->GetHistory()->GetTransform(dn - 1).TransformPoint(hitPoint) << " Down2 "
0203 << touch->GetHistory()->GetTransform(dn - 2).TransformPoint(hitPoint) << " Down3 "
0204 << touch->GetHistory()->GetTransform(dn - 3).TransformPoint(hitPoint) << " Local "
0205 << local;
0206 #endif
0207 if (fiducialCut_) {
0208 int layertype = hgcons_->layerType(layer);
0209 int frontBack = HGCalTypes::layerFrontBack(layertype);
0210 if (guardRing_->exclude(local, iz, frontBack, layer, uv.first, uv.second) ||
0211 guardRingPartial_->exclude(local, iz, frontBack, layer, uv.first, uv.second)) {
0212 id = 0;
0213 #ifdef EDM_ML_DEBUG
0214 edm::LogVerbatim("HGCSim") << "Rejected by GuardRing cutoff *****";
0215 #endif
0216 }
0217 }
0218 if ((rejectMB_) && (mouseBite_->exclude(local, iz, layer, uv.first, uv.second))) {
0219 id = 0;
0220 #ifdef EDM_ML_DEBUG
0221 edm::LogVerbatim("HGCSim") << "Rejected by MouseBite cutoff *****";
0222 #endif
0223 }
0224 }
0225 #ifdef EDM_ML_DEBUG
0226 if (id != 0)
0227 edm::LogVerbatim("HGCSim") << HGCSiliconDetId(id);
0228 #endif
0229 if ((id != 0) && checkID_) {
0230 HGCSiliconDetId hid1(id);
0231 bool cshift = (hgcons_->cassetteShiftSilicon(hid1.zside(), hid1.layer(), hid1.waferU(), hid1.waferV()));
0232 std::string_view pid = (cshift ? "HGCSim" : "HGCalSim");
0233 bool debug = (verbose_ > 0) ? true : false;
0234 auto xy = hgcons_->locateCell(hid1, debug);
0235 double xx = (hid1.zside() > 0) ? xy.first : -xy.first;
0236 double dx = xx - (hitPoint.x() / CLHEP::cm);
0237 double dy = xy.second - (hitPoint.y() / CLHEP::cm);
0238 double diff = (dx * dx + dy * dy);
0239 constexpr double tol = 2.0 * 2.0;
0240 bool valid1 = hgcons_->isValidHex8(hid1.layer(), hid1.waferU(), hid1.waferV(), hid1.cellU(), hid1.cellV(), true);
0241 if ((diff > tol) || (!valid1))
0242 pid = "HGCalError";
0243 auto partn = hgcons_->waferTypeRotation(hid1.layer(), hid1.waferU(), hid1.waferV(), false, false);
0244 int indx = HGCalWaferIndex::waferIndex(layer, hid1.waferU(), hid1.waferV());
0245 edm::LogVerbatim(pid) << "CheckID " << HGCSiliconDetId(id) << " Layer:Module:Cell:ModuleLev " << layer << ":"
0246 << module << ":" << cell << ":" << moduleLev << " SimWt:history " << useSimWt_ << ":"
0247 << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_ << " input position: ("
0248 << hitPoint.x() / CLHEP::cm << ", " << hitPoint.y() / CLHEP::cm << ":"
0249 << convertRadToDeg(std::atan2(hitPoint.y(), hitPoint.x())) << "); position from ID (" << xx
0250 << ", " << xy.second << ") distance " << dx << ":" << dy << ":" << std::sqrt(diff)
0251 << " Valid " << valid1 << " Wafer type|rotation " << partn.first << ":" << partn.second
0252 << " Part:Orient:Cassette " << std::get<1>(hgcons_->waferFileInfo(indx)) << ":"
0253 << std::get<2>(hgcons_->waferFileInfo(indx)) << ":"
0254 << std::get<3>(hgcons_->waferFileInfo(indx)) << " CassetteShift " << cshift;
0255 if ((diff > tol) || (!valid1)) {
0256 printDetectorLevels(touch);
0257 hgcons_->locateCell(hid1, true);
0258 }
0259 }
0260
0261 if ((id != 0) && calibCells_)
0262 calibCell_ = calibCell(id);
0263
0264 return id;
0265 }
0266
0267 void HGCalSD::update(const BeginOfJob* job) {
0268 if (hgcons_ != nullptr) {
0269 geom_mode_ = hgcons_->geomMode();
0270 slopeMin_ = hgcons_->minSlope();
0271 levelT1_ = hgcons_->levelTop(0);
0272 levelT2_ = hgcons_->levelTop(1);
0273 if (dd4hep_) {
0274 ++levelT1_;
0275 ++levelT2_;
0276 }
0277 useSimWt_ = hgcons_->getParameter()->useSimWt_;
0278 int useOffset = hgcons_->getParameter()->useOffset_;
0279 waferSize_ = hgcons_->waferSize(false);
0280 double mouseBite = hgcons_->mouseBite(false);
0281 guardRingOffset_ = hgcons_->guardRingOffset(false);
0282 double sensorSizeOffset = hgcons_->sensorSizeOffset(false);
0283 if (useOffset > 0) {
0284 rejectMB_ = true;
0285 fiducialCut_ = true;
0286 }
0287 double mouseBiteNew = (fiducialCut_) ? (mouseBite + guardRingOffset_ + sensorSizeOffset / cos30deg_) : mouseBite;
0288 mouseBiteCut_ = waferSize_ * tan30deg_ - mouseBiteNew;
0289 #ifdef EDM_ML_DEBUG
0290 edm::LogVerbatim("HGCSim") << "HGCalSD::Initialized with mode " << geom_mode_ << " Slope cut " << slopeMin_
0291 << " top Level " << levelT1_ << ":" << levelT2_ << " useSimWt " << useSimWt_ << " wafer "
0292 << waferSize_ << ":" << mouseBite << ":" << guardRingOffset_ << ":" << sensorSizeOffset
0293 << ":" << mouseBiteNew << ":" << mouseBiteCut_ << " useOffset " << useOffset
0294 << " dd4hep " << dd4hep_;
0295 #endif
0296
0297 numberingScheme_ = std::make_unique<HGCalNumberingScheme>(*hgcons_, mydet_, nameX_, missingFile_);
0298 if (rejectMB_)
0299 mouseBite_ = std::make_unique<HGCMouseBite>(*hgcons_, angles_, mouseBiteCut_, waferRot_);
0300 if (fiducialCut_) {
0301 guardRing_ = std::make_unique<HGCGuardRing>(*hgcons_);
0302 guardRingPartial_ = std::make_unique<HGCGuardRingPartial>(*hgcons_);
0303 }
0304
0305
0306 calibCellRHD_ = hgcons_->calibCellRad(true);
0307 calibCellRLD_ = hgcons_->calibCellRad(false);
0308 calibCellFullHD_ = hgcons_->calibCells(true, true);
0309 calibCellPartHD_ = hgcons_->calibCells(true, false);
0310 calibCellFullLD_ = hgcons_->calibCells(false, true);
0311 calibCellPartLD_ = hgcons_->calibCells(false, false);
0312 calibCells_ = ((!calibCellFullHD_.empty()) || (!calibCellPartHD_.empty()));
0313 calibCells_ |= ((!calibCellFullLD_.empty()) || (!calibCellPartLD_.empty()));
0314 #ifdef EDM_ML_DEBUG
0315 edm::LogVerbatim("HGCSim") << "HGCalSD::Calibration cells initialized with flag " << calibCells_;
0316 edm::LogVerbatim("HGCSim") << " Radii " << calibCellRHD_ << ":" << calibCellRLD_;
0317 std::ostringstream st1;
0318 for (const auto& cell : calibCellFullHD_)
0319 st1 << " " << cell;
0320 edm::LogVerbatim("HGCSim") << calibCellFullHD_.size() << " cells for High Density full wafers: " << st1.str();
0321 std::ostringstream st2;
0322 for (const auto& cell : calibCellPartHD_)
0323 st2 << " " << cell;
0324 edm::LogVerbatim("HGCSim") << calibCellPartHD_.size() << " cells for High Density partial wafers: " << st2.str();
0325 std::ostringstream st3;
0326 for (const auto& cell : calibCellFullLD_)
0327 st3 << " " << cell;
0328 edm::LogVerbatim("HGCSim") << calibCellFullLD_.size() << " cells for Low Density full wafers: " << st3.str();
0329 std::ostringstream st4;
0330 for (const auto& cell : calibCellPartLD_)
0331 st4 << " " << cell;
0332 edm::LogVerbatim("HGCSim") << calibCellPartLD_.size() << " cells for Low Density partial wafers: " << st4.str();
0333 #endif
0334 } else {
0335 throw cms::Exception("Unknown", "HGCalSD") << "Cannot find HGCalDDDConstants for " << nameX_ << "\n";
0336 }
0337 if ((nHC_ > 1) && calibCells_) {
0338 newCollection(collName_[1], ps_);
0339 cellOffset_ = std::make_unique<HGCalCellOffset>(
0340 waferSize_, hgcons_->getUVMax(0), hgcons_->getUVMax(1), guardRingOffset_, mouseBiteCut_);
0341 }
0342 }
0343
0344 void HGCalSD::initRun() {}
0345
0346 bool HGCalSD::filterHit(CaloG4Hit* aHit, double time) {
0347 return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit_));
0348 }
0349
0350 void HGCalSD::processSecondHit(const G4Step* aStep, const G4Track* theTrack) {
0351 if (calibCell_) {
0352 float edEM(edepositEM), edHad(edepositHAD);
0353 currentID[1] = currentID[0];
0354 edepositEM *= fraction_;
0355 edepositHAD *= fraction_;
0356 if (!hitExists(aStep, 1)) {
0357 currentHit[1] = createNewHit(aStep, theTrack, 1);
0358 }
0359 edepositEM = edEM;
0360 edepositHAD = edHad;
0361 }
0362 }
0363
0364 uint32_t HGCalSD::setDetUnitId(int layer, int module, int cell, int iz, G4ThreeVector& pos) {
0365 uint32_t id = numberingScheme_ ? numberingScheme_->getUnitID(layer, module, cell, iz, pos, weight_) : 0;
0366 if (cornerMinMask_ > 2) {
0367 if (hgcons_->maskCell(DetId(id), cornerMinMask_)) {
0368 id = 0;
0369 ignoreRejection();
0370 }
0371 }
0372 if (hgcons_->waferHexagon8File() || (id == 0))
0373 ignoreRejection();
0374
0375 return id;
0376 }
0377
0378 bool HGCalSD::calibCell(const uint32_t& id) {
0379 bool flag(false);
0380 int type, zside, layer, waferU, waferV, cellU, cellV;
0381 HGCSiliconDetId(id).unpack(type, zside, layer, waferU, waferV, cellU, cellV);
0382 HGCalParameters::waferInfo info = hgcons_->waferInfo(layer, waferU, waferV);
0383 bool hd = HGCalTypes::waferHD(info.type);
0384 bool full = HGCalTypes::waferFull(info.part);
0385 int indx = 100 * cellU + cellV;
0386 if (hd) {
0387 if (full)
0388 flag = (std::find(calibCellFullHD_.begin(), calibCellFullHD_.end(), indx) != calibCellFullHD_.end());
0389 else
0390 flag = (std::find(calibCellPartHD_.begin(), calibCellPartHD_.end(), indx) != calibCellPartHD_.end());
0391 } else {
0392 if (full)
0393 flag = (std::find(calibCellFullLD_.begin(), calibCellFullLD_.end(), indx) != calibCellFullLD_.end());
0394 else
0395 flag = (std::find(calibCellPartLD_.begin(), calibCellPartLD_.end(), indx) != calibCellPartLD_.end());
0396 }
0397 if (flag) {
0398 int32_t place =
0399 HGCalCell::cellPlacementIndex(zside, HGCalTypes::layerFrontBack(hgcons_->layerType(layer)), info.orient);
0400 int32_t type = hd ? 0 : 1;
0401 double num = hd ? (M_PI * calibCellRHD_ * calibCellRHD_) : (M_PI * calibCellRLD_ * calibCellRLD_);
0402 double bot = cellOffset_->cellAreaUV(cellU, cellV, place, type, true);
0403 fraction_ = (bot > 0 && bot > num) ? (num / bot) : 1.0;
0404 #ifdef EDM_ML_DEBUG
0405 edm::LogVerbatim("HGCSim") << HGCSiliconDetId(id) << " CalibrationCell flag " << flag << " and fraction " << num
0406 << ":" << bot << ":" << fraction_;
0407 #endif
0408 } else {
0409 #ifdef EDM_ML_DEBUG
0410 edm::LogVerbatim("HGCSim") << HGCSiliconDetId(id) << " CalibrationCell flag " << flag;
0411 #endif
0412 }
0413 return flag;
0414 }