File indexing completed on 2025-01-31 09:21:47
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 #include "G4HCofThisEvent.hh"
0044 #include "G4Step.hh"
0045 #include "G4ThreeVector.hh"
0046 #include "G4SDManager.hh"
0047 #include "G4ios.hh"
0048 #include "G4Track.hh"
0049 #ifdef WITH_G4OPTICKS
0050 # include "G4Opticks.hh"
0051 # include "TrackInfo.hh"
0052 # include "OpticksGenstep.h"
0053 # include "OpticksFlags.hh"
0054 # include "G4OpticksHit.hh"
0055 # include "G4Cerenkov.hh"
0056 # include "G4Event.hh"
0057 # include "G4MaterialPropertiesTable.hh"
0058 # include "G4PhysicalConstants.hh"
0059 # include "G4RunManager.hh"
0060 # include "G4SteppingManager.hh"
0061 # include "G4SystemOfUnits.hh"
0062 # include "G4UnitsTable.hh"
0063 # include "G4VProcess.hh"
0064 # include "G4VRestDiscreteProcess.hh"
0065 # include "PhotonSD.hh"
0066 # include "G4Cerenkov.hh"
0067 # include "G4Scintillation.hh"
0068 # include "G4Version.hh"
0069 #endif
0070
0071 #include "lArTPCSD.hh"
0072 #include "ConfigurationManager.hh"
0073
0074 lArTPCSD::lArTPCSD(G4String name)
0075 : G4VSensitiveDetector(name)
0076 {
0077 G4String HCname = name + "_HC";
0078 collectionName.insert(HCname);
0079 verbose = ConfigurationManager::getInstance()->isEnable_verbose();
0080 if(verbose)
0081 {
0082 G4cout << collectionName.size() << " lArTPCSD name: " << name
0083 << " collection Name: " << HCname << G4endl;
0084 }
0085 fHCID = -1;
0086 first = true;
0087 }
0088
0089 void lArTPCSD::Initialize(G4HCofThisEvent* hce)
0090 {
0091 flArTPCHitsCollection =
0092 new lArTPCHitsCollection(SensitiveDetectorName, collectionName[0]);
0093 if(fHCID < 0)
0094 {
0095 if(verbose)
0096 {
0097 G4cout << "lArTPCSD::Initialize: " << SensitiveDetectorName << " "
0098 << collectionName[0] << G4endl;
0099 }
0100 fHCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
0101 }
0102 hce->AddHitsCollection(fHCID, flArTPCHitsCollection);
0103 }
0104
0105 G4bool lArTPCSD::ProcessHits(G4Step* aStep, G4TouchableHistory*)
0106 {
0107 G4double edep = aStep->GetTotalEnergyDeposit();
0108 if(edep == 0.)
0109 return false;
0110
0111 G4Track* aTrack = aStep->GetTrack();
0112 G4double charge = aTrack->GetDynamicParticle()->GetCharge();
0113 if(charge == 0)
0114 return false;
0115 G4double ds = aStep->GetStepLength();
0116 lArTPCHit* newHit = new lArTPCHit(
0117 NumElectrons(edep, ds), aStep->GetPostStepPoint()->GetPosition().getX(),
0118 aStep->GetPostStepPoint()->GetPosition().getY(),
0119 aStep->GetPostStepPoint()->GetPosition().getZ());
0120 flArTPCHitsCollection->insert(newHit);
0121 #ifdef WITH_G4OPTICKS
0122 if(ConfigurationManager::getInstance()->isEnable_opticks())
0123 {
0124 if(first)
0125 {
0126 aMaterial = aTrack->GetMaterial();
0127 materialIndex = aMaterial->GetIndex();
0128 if(verbose)
0129 {
0130 G4cout << "*******************************" << G4endl;
0131 G4cout << "RadiatorSD::ProcessHits initializing Material: "
0132 << aMaterial->GetName() << " " << G4endl;
0133 G4cout << "RadiatorSD::ProcessHits: Name "
0134 << aStep->GetPreStepPoint()
0135 ->GetPhysicalVolume()
0136 ->GetLogicalVolume()
0137 ->GetName()
0138 << G4endl;
0139 }
0140 aMaterialPropertiesTable = aMaterial->GetMaterialPropertiesTable();
0141 if(verbose)
0142 {
0143 aMaterialPropertiesTable->DumpTable();
0144 }
0145
0146
0147
0148 # if(G4VERSION_NUMBER > 1072)
0149 YieldRatio =
0150 aMaterialPropertiesTable->GetConstProperty(kSCINTILLATIONYIELD1) /
0151 aMaterialPropertiesTable->GetConstProperty(
0152 kSCINTILLATIONYIELD2);
0153 FastTimeConstant = aMaterialPropertiesTable->GetConstProperty(
0154 kSCINTILLATIONTIMECONSTANT1);
0155 SlowTimeConstant = aMaterialPropertiesTable->GetConstProperty(
0156 kSCINTILLATIONTIMECONSTANT2);
0157 # else
0158 Fast_Intensity = aMaterialPropertiesTable->GetProperty(kFASTCOMPONENT);
0159 Slow_Intensity = aMaterialPropertiesTable->GetProperty(kSLOWCOMPONENT);
0160 YieldRatio = aMaterialPropertiesTable->GetConstProperty(kYIELDRATIO);
0161 # endif
0162 ScintillationType = Slow;
0163
0164
0165
0166 Rindex = aMaterialPropertiesTable->GetProperty("RINDEX");
0167 # if(G4VERSION_NUMBER > 1072)
0168 Pmin = Rindex->GetMinEnergy();
0169 Pmax = Rindex->GetMaxEnergy();
0170 # else
0171 Pmin = Rindex->GetMinLowEdgeEnergy();
0172 Pmax = Rindex->GetMaxLowEdgeEnergy();
0173 # endif
0174 dp = Pmax - Pmin;
0175 if(verbose)
0176 {
0177 G4cout << "nMax: " << nMax << "Pmin: " << Pmin << "Pmax: " << Pmax
0178 << "dp: " << dp << G4endl;
0179 Rindex->DumpValues();
0180 }
0181
0182 first = false;
0183 }
0184 G4int Sphotons = 0;
0185 G4int Cphotons = 0;
0186
0187
0188
0189 G4double maxCos = 0.0;
0190 G4double maxSin2 = 0.0;
0191 G4double beta = 0.0;
0192 G4double beta1 = 0.0;
0193 G4double beta2 = 0.0;
0194 G4double BetaInverse = 0.0;
0195 G4double MeanNumberOfPhotons1 = 0.0;
0196 G4double MeanNumberOfPhotons2 = 0.0;
0197 G4SteppingManager* fpSteppingManager = G4EventManager::GetEventManager()
0198 ->GetTrackingManager()
0199 ->GetSteppingManager();
0200 G4StepStatus stepStatus = fpSteppingManager->GetfStepStatus();
0201 if(stepStatus != fAtRestDoItProc)
0202 {
0203 G4ProcessVector* procPost = fpSteppingManager->GetfPostStepDoItVector();
0204 size_t MAXofPostStepLoops = fpSteppingManager->GetMAXofPostStepLoops();
0205 for(size_t i3 = 0; i3 < MAXofPostStepLoops; i3++)
0206 {
0207 if((*procPost)[i3]->GetProcessName() == "Cerenkov")
0208 {
0209 G4Cerenkov* proc = (G4Cerenkov*) (*procPost)[i3];
0210 thePhysicsTable = proc->GetPhysicsTable();
0211 CerenkovAngleIntegrals =
0212 (G4PhysicsOrderedFreeVector*) ((*thePhysicsTable)(materialIndex));
0213 Cphotons = proc->GetNumPhotons();
0214 if(Cphotons > 0)
0215 {
0216 beta1 = aStep->GetPreStepPoint()->GetBeta();
0217 beta2 = aStep->GetPostStepPoint()->GetBeta();
0218 beta = (beta1 + beta2) * 0.5;
0219 BetaInverse = 1. / beta;
0220 maxCos = BetaInverse / nMax;
0221 maxSin2 = (1.0 - maxCos) * (1.0 + maxCos);
0222 MeanNumberOfPhotons1 =
0223 proc->GetAverageNumberOfPhotons(charge, beta1, aMaterial, Rindex);
0224 MeanNumberOfPhotons2 =
0225 proc->GetAverageNumberOfPhotons(charge, beta2, aMaterial, Rindex);
0226 }
0227 }
0228 if((*procPost)[i3]->GetProcessName() == "Scintillation")
0229 {
0230 G4Scintillation* proc1 = (G4Scintillation*) (*procPost)[i3];
0231 Sphotons = proc1->GetNumPhotons();
0232 }
0233 }
0234 }
0235 tSphotons += Sphotons;
0236 tCphotons += Cphotons;
0237 G4ThreeVector deltaPosition = aStep->GetDeltaPosition();
0238 G4double ScintillationTime = 0. * ns;
0239 G4int scntId = 1;
0240 G4StepPoint* pPreStepPoint = aStep->GetPreStepPoint();
0241 G4ThreeVector x0 = pPreStepPoint->GetPosition();
0242 G4ThreeVector p0 = aStep->GetDeltaPosition().unit();
0243
0244
0245
0246 if(Sphotons > 0)
0247 {
0248 G4double ScintillationRiseTime = 0.0;
0249 G4Opticks::Get()->collectGenstep_G4Scintillation_1042(
0250 aTrack, aStep, Sphotons, scntId, ScintillationTime,
0251 ScintillationRiseTime);
0252 }
0253
0254
0255
0256 if(Cphotons > 0)
0257 {
0258 G4Opticks::Get()->collectGenstep_G4Cerenkov_1042(
0259 aTrack, aStep, Cphotons, BetaInverse, Pmin, Pmax, maxCos, maxSin2,
0260 MeanNumberOfPhotons1, MeanNumberOfPhotons2);
0261 }
0262 G4Opticks* g4ok = G4Opticks::Get();
0263 G4RunManager* rm = G4RunManager::GetRunManager();
0264 const G4Event* event = rm->GetCurrentEvent();
0265 G4int eventid = event->GetEventID();
0266 G4OpticksHit hit;
0267 unsigned num_photons = g4ok->getNumPhotons();
0268 if(num_photons > ConfigurationManager::getInstance()->getMaxPhotons())
0269 {
0270 g4ok->propagateOpticalPhotons(eventid);
0271 G4HCtable* hctable = G4SDManager::GetSDMpointer()->GetHCtable();
0272 for(G4int i = 0; i < hctable->entries(); ++i)
0273 {
0274 std::string sdn = hctable->GetSDname(i);
0275 std::size_t found = sdn.find("Photondetector");
0276 if(found != std::string::npos)
0277 {
0278 PhotonSD* aSD =
0279 (PhotonSD*) G4SDManager::GetSDMpointer()->FindSensitiveDetector(
0280 sdn);
0281 aSD->AddOpticksHits();
0282 }
0283 }
0284 g4ok->reset();
0285 }
0286 }
0287 #endif
0288 return true;
0289 }
0290
0291 void lArTPCSD::EndOfEvent(G4HCofThisEvent*)
0292 {
0293 tSphotons = 0;
0294 tCphotons = 0;
0295 G4int NbHits = flArTPCHitsCollection->entries();
0296 if(verbose)
0297 {
0298 G4cout << " Number of lArTPCHits: " << NbHits << G4endl;
0299 }
0300 }
0301
0302 G4double lArTPCSD::NumElectrons(G4double edep, G4double ds)
0303 {
0304
0305 G4double fGeVToElectrons = 4.237e+07;
0306 G4double fModBoxA = 0.930;
0307 G4double fModBoxB = 0.212;
0308 G4double EFieldStep = 0.5;
0309 G4double recomb = 0.0;
0310 G4double dEdx = (ds <= 0.0) ? 0.0 : edep / ds;
0311 if(dEdx < 1.)
0312 dEdx = 1.;
0313 if(ds > 0)
0314 {
0315 G4double Xi = fModBoxB * dEdx / EFieldStep;
0316 recomb = std::log(fModBoxA + Xi) / Xi;
0317 }
0318 else
0319 {
0320 recomb = 0.0;
0321 }
0322 G4double fNumIonElectrons = fGeVToElectrons * 1.e-3 * edep * recomb;
0323 return fNumIonElectrons;
0324 }