File indexing completed on 2024-06-13 03:24:12
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 #ifdef EDM_ML_DEBUG
0194 edm::LogVerbatim("HGCSim") << "ID Layer " << layer << " Module " << module << " Cell " << cell << " " << std::hex
0195 << id << std::dec << " " << HGCSiliconDetId(id);
0196 #endif
0197 if ((rejectMB_ || fiducialCut_) && id != 0) {
0198 auto uv = HGCSiliconDetId(id).waferUV();
0199 #ifdef EDM_ML_DEBUG
0200 edm::LogVerbatim("HGCSim") << "ID " << std::hex << id << std::dec << " " << HGCSiliconDetId(id);
0201 #endif
0202 G4ThreeVector local = (touch->GetHistory()->GetTransform(dn - moduleLev).TransformPoint(hitPoint));
0203 #ifdef EDM_ML_DEBUG
0204 edm::LogVerbatim("HGCSim") << "Global Point " << hitPoint << " Down0 "
0205 << touch->GetHistory()->GetTransform(dn).TransformPoint(hitPoint) << " Down1 "
0206 << touch->GetHistory()->GetTransform(dn - 1).TransformPoint(hitPoint) << " Down2 "
0207 << touch->GetHistory()->GetTransform(dn - 2).TransformPoint(hitPoint) << " Down3 "
0208 << touch->GetHistory()->GetTransform(dn - 3).TransformPoint(hitPoint) << " Local "
0209 << local;
0210 #endif
0211 if (fiducialCut_) {
0212 int layertype = hgcons_->layerType(layer);
0213 int frontBack = HGCalTypes::layerFrontBack(layertype);
0214 if (guardRing_->exclude(local, iz, frontBack, layer, uv.first, uv.second) ||
0215 guardRingPartial_->exclude(local, iz, frontBack, layer, uv.first, uv.second)) {
0216 id = 0;
0217 #ifdef EDM_ML_DEBUG
0218 edm::LogVerbatim("HGCSim") << "Rejected by GuardRing cutoff *****";
0219 #endif
0220 }
0221 }
0222 if ((rejectMB_) && (mouseBite_->exclude(local, iz, layer, uv.first, uv.second))) {
0223 id = 0;
0224 #ifdef EDM_ML_DEBUG
0225 edm::LogVerbatim("HGCSim") << "Rejected by MouseBite cutoff *****";
0226 #endif
0227 }
0228 }
0229 #ifdef EDM_ML_DEBUG
0230 if (id != 0)
0231 edm::LogVerbatim("HGCSim") << HGCSiliconDetId(id);
0232 #endif
0233 if ((id != 0) && checkID_) {
0234 HGCSiliconDetId hid1(id);
0235 bool cshift = (hgcons_->cassetteShiftSilicon(hid1.zside(), hid1.layer(), hid1.waferU(), hid1.waferV()));
0236 std::string_view pid = (cshift ? "HGCSim" : "HGCalSim");
0237 bool debug = (verbose_ > 0) ? true : false;
0238 auto xy = hgcons_->locateCell(hid1, debug);
0239 double xx = (hid1.zside() > 0) ? xy.first : -xy.first;
0240 double dx = xx - (hitPoint.x() / CLHEP::cm);
0241 double dy = xy.second - (hitPoint.y() / CLHEP::cm);
0242 double diff = (dx * dx + dy * dy);
0243 constexpr double tol = 2.0 * 2.0;
0244 bool valid1 = hgcons_->isValidHex8(hid1.layer(), hid1.waferU(), hid1.waferV(), hid1.cellU(), hid1.cellV(), true);
0245 if ((diff > tol) || (!valid1))
0246 pid = "HGCalError";
0247 auto partn = hgcons_->waferTypeRotation(hid1.layer(), hid1.waferU(), hid1.waferV(), false, false);
0248 int indx = HGCalWaferIndex::waferIndex(layer, hid1.waferU(), hid1.waferV());
0249 edm::LogVerbatim(pid) << "CheckID " << HGCSiliconDetId(id) << " Layer:Module:Cell:ModuleLev " << layer << ":"
0250 << module << ":" << cell << ":" << moduleLev << " SimWt:history " << useSimWt_ << ":"
0251 << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_ << " input position: ("
0252 << hitPoint.x() / CLHEP::cm << ", " << hitPoint.y() / CLHEP::cm << ":"
0253 << convertRadToDeg(std::atan2(hitPoint.y(), hitPoint.x())) << "); position from ID (" << xx
0254 << ", " << xy.second << ") distance " << dx << ":" << dy << ":" << std::sqrt(diff)
0255 << " Valid " << valid1 << " Wafer type|rotation " << partn.first << ":" << partn.second
0256 << " Part:Orient:Cassette " << std::get<1>(hgcons_->waferFileInfo(indx)) << ":"
0257 << std::get<2>(hgcons_->waferFileInfo(indx)) << ":"
0258 << std::get<3>(hgcons_->waferFileInfo(indx)) << " CassetteShift " << cshift;
0259 if ((diff > tol) || (!valid1)) {
0260 printDetectorLevels(touch);
0261 hgcons_->locateCell(hid1, true);
0262 }
0263 }
0264
0265 if ((id != 0) && calibCells_)
0266 calibCell_ = calibCell(id);
0267
0268 return id;
0269 }
0270
0271 void HGCalSD::update(const BeginOfJob* job) {
0272 if (hgcons_ != nullptr) {
0273 geom_mode_ = hgcons_->geomMode();
0274 slopeMin_ = hgcons_->minSlope();
0275 levelT1_ = hgcons_->levelTop(0);
0276 levelT2_ = hgcons_->levelTop(1);
0277 if (dd4hep_) {
0278 ++levelT1_;
0279 ++levelT2_;
0280 }
0281 useSimWt_ = hgcons_->getParameter()->useSimWt_;
0282 int useOffset = hgcons_->getParameter()->useOffset_;
0283 waferSize_ = hgcons_->waferSize(false);
0284 double mouseBite = hgcons_->mouseBite(false);
0285 guardRingOffset_ = hgcons_->guardRingOffset(false);
0286 sensorSizeOffset_ = hgcons_->sensorSizeOffset(false);
0287 if (useOffset > 0) {
0288 rejectMB_ = true;
0289 fiducialCut_ = true;
0290 }
0291 double mouseBiteNew = (fiducialCut_) ? (mouseBite + guardRingOffset_ + sensorSizeOffset_ / cos30deg_) : mouseBite;
0292 mouseBiteCut_ = waferSize_ * tan30deg_ - mouseBiteNew;
0293 #ifdef EDM_ML_DEBUG
0294 edm::LogVerbatim("HGCSim") << "HGCalSD::Initialized with mode " << geom_mode_ << " Slope cut " << slopeMin_
0295 << " top Level " << levelT1_ << ":" << levelT2_ << " useSimWt " << useSimWt_ << " wafer "
0296 << waferSize_ << ":" << mouseBite << ":" << guardRingOffset_ << ":" << sensorSizeOffset_
0297 << ":" << mouseBiteNew << ":" << mouseBiteCut_ << " useOffset " << useOffset
0298 << " dd4hep " << dd4hep_;
0299 #endif
0300
0301 numberingScheme_ = std::make_unique<HGCalNumberingScheme>(*hgcons_, mydet_, nameX_, missingFile_);
0302 if (rejectMB_)
0303 mouseBite_ = std::make_unique<HGCMouseBite>(*hgcons_, angles_, mouseBiteCut_, waferRot_);
0304 if (fiducialCut_) {
0305 guardRing_ = std::make_unique<HGCGuardRing>(*hgcons_);
0306 guardRingPartial_ = std::make_unique<HGCGuardRingPartial>(*hgcons_);
0307 }
0308
0309
0310 calibCellRHD_ = hgcons_->calibCellRad(true);
0311 calibCellRLD_ = hgcons_->calibCellRad(false);
0312 calibCellFullHD_ = hgcons_->calibCells(true, true);
0313 calibCellPartHD_ = hgcons_->calibCells(true, false);
0314 calibCellFullLD_ = hgcons_->calibCells(false, true);
0315 calibCellPartLD_ = hgcons_->calibCells(false, false);
0316 calibCells_ = ((!calibCellFullHD_.empty()) || (!calibCellPartHD_.empty()));
0317 calibCells_ |= ((!calibCellFullLD_.empty()) || (!calibCellPartLD_.empty()));
0318 #ifdef EDM_ML_DEBUG
0319 edm::LogVerbatim("HGCSim") << "HGCalSD::Calibration cells initialized with flag " << calibCells_;
0320 edm::LogVerbatim("HGCSim") << " Radii " << calibCellRHD_ << ":" << calibCellRLD_;
0321 std::ostringstream st1;
0322 for (const auto& cell : calibCellFullHD_)
0323 st1 << " " << cell;
0324 edm::LogVerbatim("HGCSim") << calibCellFullHD_.size() << " cells for High Density full wafers: " << st1.str();
0325 std::ostringstream st2;
0326 for (const auto& cell : calibCellPartHD_)
0327 st2 << " " << cell;
0328 edm::LogVerbatim("HGCSim") << calibCellPartHD_.size() << " cells for High Density partial wafers: " << st2.str();
0329 std::ostringstream st3;
0330 for (const auto& cell : calibCellFullLD_)
0331 st3 << " " << cell;
0332 edm::LogVerbatim("HGCSim") << calibCellFullLD_.size() << " cells for Low Density full wafers: " << st3.str();
0333 std::ostringstream st4;
0334 for (const auto& cell : calibCellPartLD_)
0335 st4 << " " << cell;
0336 edm::LogVerbatim("HGCSim") << calibCellPartLD_.size() << " cells for Low Density partial wafers: " << st4.str();
0337 #endif
0338 } else {
0339 throw cms::Exception("Unknown", "HGCalSD") << "Cannot find HGCalDDDConstants for " << nameX_ << "\n";
0340 }
0341 if ((nHC_ > 1) && calibCells_)
0342 newCollection(collName_[1], ps_);
0343 cellOffset_ = std::make_unique<HGCalCellOffset>(
0344 waferSize_, hgcons_->getUVMax(0), hgcons_->getUVMax(1), guardRingOffset_, mouseBiteCut_, sensorSizeOffset_);
0345 }
0346
0347 void HGCalSD::initRun() {}
0348
0349 bool HGCalSD::filterHit(CaloG4Hit* aHit, double time) {
0350 return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit_));
0351 }
0352
0353 void HGCalSD::processSecondHit(const G4Step* aStep, const G4Track* theTrack) {
0354 if (calibCell_) {
0355 float edEM(edepositEM), edHad(edepositHAD);
0356 currentID[1] = currentID[0];
0357 edepositEM *= fraction_;
0358 edepositHAD *= fraction_;
0359 if (!hitExists(aStep, 1)) {
0360 currentHit[1] = createNewHit(aStep, theTrack, 1);
0361 }
0362 edepositEM = edEM;
0363 edepositHAD = edHad;
0364 }
0365 }
0366
0367 uint32_t HGCalSD::setDetUnitId(int layer, int module, int cell, int iz, G4ThreeVector& pos) {
0368 uint32_t id = numberingScheme_ ? numberingScheme_->getUnitID(layer, module, cell, iz, pos, weight_) : 0;
0369 if (cornerMinMask_ > 2) {
0370 if (hgcons_->maskCell(DetId(id), cornerMinMask_)) {
0371 id = 0;
0372 ignoreRejection();
0373 }
0374 }
0375 if (hgcons_->waferHexagon8File() || (id == 0))
0376 ignoreRejection();
0377
0378 return id;
0379 }
0380
0381 bool HGCalSD::calibCell(const uint32_t& id) {
0382 bool flag(false);
0383 int type, zside, layer, waferU, waferV, cellU, cellV;
0384 HGCSiliconDetId(id).unpack(type, zside, layer, waferU, waferV, cellU, cellV);
0385 HGCalParameters::waferInfo info = hgcons_->waferInfo(layer, waferU, waferV);
0386 bool hd = HGCalTypes::waferHD(info.type);
0387 bool full = HGCalTypes::waferFull(info.part);
0388 int indx = 100 * cellU + cellV;
0389 if (hd) {
0390 if (full)
0391 flag = (std::find(calibCellFullHD_.begin(), calibCellFullHD_.end(), indx) != calibCellFullHD_.end());
0392 else
0393 flag = (std::find(calibCellPartHD_.begin(), calibCellPartHD_.end(), indx) != calibCellPartHD_.end());
0394 } else {
0395 if (full)
0396 flag = (std::find(calibCellFullLD_.begin(), calibCellFullLD_.end(), indx) != calibCellFullLD_.end());
0397 else
0398 flag = (std::find(calibCellPartLD_.begin(), calibCellPartLD_.end(), indx) != calibCellPartLD_.end());
0399 }
0400 if (flag) {
0401 int32_t place =
0402 HGCalCell::cellPlacementIndex(zside, HGCalTypes::layerFrontBack(hgcons_->layerType(layer)), info.orient);
0403 int32_t type = hd ? 0 : 1;
0404 double num = hd ? (M_PI * calibCellRHD_ * calibCellRHD_) : (M_PI * calibCellRLD_ * calibCellRLD_);
0405 double bot = cellOffset_->cellAreaUV(cellU, cellV, place, type, true);
0406 fraction_ = (bot > 0 && bot > num) ? (num / bot) : 1.0;
0407 #ifdef EDM_ML_DEBUG
0408 edm::LogVerbatim("HGCSim") << HGCSiliconDetId(id) << " CalibrationCell flag " << flag << " and fraction " << num
0409 << ":" << bot << ":" << fraction_;
0410 #endif
0411 } else {
0412 #ifdef EDM_ML_DEBUG
0413 edm::LogVerbatim("HGCSim") << HGCSiliconDetId(id) << " CalibrationCell flag " << flag;
0414 #endif
0415 }
0416 return flag;
0417 }