Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 // 27 /// \file Par02FastSimModelHCal.cc 28 /// \brief Implementation of the Par02FastSimModelHCal class 29 30 #include "Par02FastSimModelHCal.hh" 31 32 #include "Par02EventInformation.hh" 33 #include "Par02Output.hh" 34 #include "Par02PrimaryParticleInformation.hh" 35 #include "Par02Smearer.hh" 36 37 #include "G4AnalysisManager.hh" 38 #include "G4Event.hh" 39 #include "G4RunManager.hh" 40 #include "G4SystemOfUnits.hh" 41 #include "G4Track.hh" 42 #include "Randomize.hh" 43 44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 45 46 Par02FastSimModelHCal::Par02FastSimModelHCal(G4String aModelName, G4Region* aEnvelope, 47 Par02DetectorParametrisation::Parametrisation aType) 48 : G4VFastSimulationModel(aModelName, aEnvelope), 49 fCalculateParametrisation(), 50 fParametrisation(aType) 51 {} 52 53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 54 55 Par02FastSimModelHCal::Par02FastSimModelHCal(G4String aModelName, G4Region* aEnvelope) 56 : G4VFastSimulationModel(aModelName, aEnvelope), 57 fCalculateParametrisation(), 58 fParametrisation(Par02DetectorParametrisation::eCMS) 59 {} 60 61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 62 63 Par02FastSimModelHCal::Par02FastSimModelHCal(G4String aModelName) 64 : G4VFastSimulationModel(aModelName), 65 fCalculateParametrisation(), 66 fParametrisation(Par02DetectorParametrisation::eCMS) 67 {} 68 69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 70 71 Par02FastSimModelHCal::~Par02FastSimModelHCal() = default; 72 73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 74 75 G4bool Par02FastSimModelHCal::IsApplicable(const G4ParticleDefinition& aParticleType) 76 { 77 G4bool isOk = false; 78 // Applicable to all hadrons, i.e. any particle made of quarks 79 if (aParticleType.GetQuarkContent(1) + aParticleType.GetQuarkContent(2) 80 + aParticleType.GetQuarkContent(3) + aParticleType.GetQuarkContent(4) 81 + aParticleType.GetQuarkContent(5) + aParticleType.GetQuarkContent(6) 82 != 0) 83 { 84 isOk = true; 85 } 86 return isOk; 87 } 88 89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 90 91 G4bool Par02FastSimModelHCal::ModelTrigger(const G4FastTrack& /*aFastTrack*/) 92 { 93 return true; // No kinematical restrictions to apply the parametrisation 94 } 95 96 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 97 98 void Par02FastSimModelHCal::DoIt(const G4FastTrack& aFastTrack, G4FastStep& aFastStep) 99 { 100 // G4cout << " ________HCal model triggered _________" << G4endl; 101 102 // Kill the parameterised particle at the entrance of the hadronic calorimeter 103 aFastStep.KillPrimaryTrack(); 104 aFastStep.ProposePrimaryTrackPathLength(0.0); 105 G4double Edep = aFastTrack.GetPrimaryTrack()->GetKineticEnergy(); 106 107 // Consider only primary tracks (do nothing else for secondary hadrons) 108 G4ThreeVector Pos = aFastTrack.GetPrimaryTrack()->GetPosition(); 109 if (!aFastTrack.GetPrimaryTrack()->GetParentID()) { 110 auto info = (Par02EventInformation*)G4EventManager::GetEventManager()->GetUserInformation(); 111 if (info->GetDoSmearing()) { 112 // Smearing according to the hadronic calorimeter resolution 113 G4ThreeVector Porg = aFastTrack.GetPrimaryTrack()->GetMomentum(); 114 G4double res = fCalculateParametrisation->GetResolution(Par02DetectorParametrisation::eHCAL, 115 fParametrisation, Porg.mag()); 116 G4double eff = fCalculateParametrisation->GetEfficiency(Par02DetectorParametrisation::eHCAL, 117 fParametrisation, Porg.mag()); 118 G4double Esm; 119 Esm = std::abs(Par02Smearer::Instance()->SmearEnergy(aFastTrack.GetPrimaryTrack(), res)); 120 Par02Output::Instance()->FillHistogram(2, (Esm / MeV) / (Edep / MeV)); 121 // Setting the values of Pos, Esm, res and eff 122 auto primaryInfo = static_cast<Par02PrimaryParticleInformation*>( 123 (aFastTrack.GetPrimaryTrack()->GetDynamicParticle()->GetPrimaryParticle()) 124 ->GetUserInformation()); 125 primaryInfo->SetHCalPosition(Pos); 126 primaryInfo->SetHCalEnergy(Esm); 127 primaryInfo->SetHCalResolution(res); 128 primaryInfo->SetHCalEfficiency(eff); 129 // The (smeared) energy of the particle is deposited in the step 130 // (which corresponds to the entrance of the hadronic calorimeter) 131 aFastStep.ProposeTotalEnergyDeposited(Esm); 132 } 133 else { 134 // No smearing: simply setting the value of Edep 135 ((Par02PrimaryParticleInformation*)(const_cast<G4PrimaryParticle*>( 136 aFastTrack.GetPrimaryTrack() 137 ->GetDynamicParticle() 138 ->GetPrimaryParticle()) 139 ->GetUserInformation())) 140 ->SetHCalEnergy(Edep); 141 // The (initial) energy of the particle is deposited in the step 142 // (which corresponds to the entrance of the hadronic calorimeter) 143 aFastStep.ProposeTotalEnergyDeposited(Edep); 144 } 145 } 146 } 147 148 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 149