Loading...
Searching...
No Matches
Par04MLFastSimModel.cc
Go to the documentation of this file.
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#ifdef USE_INFERENCE
28#include <stddef.h> // for size_t
29#include <G4FastStep.hh> // for G4FastStep
30#include <G4FastTrack.hh> // for G4FastTrack
31#include <G4Track.hh> // for G4Track
32#include <G4VFastSimulationModel.hh> // for G4VFastSimulationModel
33#include "G4Electron.hh" // for G4Electron
34#include "G4FastHit.hh" // for G4FastHit
35#include "G4FastSimHitMaker.hh" // for G4FastSimHitMaker
36#include "G4Gamma.hh" // for G4Gamma
37#include "G4Positron.hh" // for G4Positron
38#include "Par04InferenceSetup.hh" // for Par04InferenceSetup
40class G4Region;
41
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43
44Par04MLFastSimModel::Par04MLFastSimModel(G4String aModelName, G4Region* aEnvelope)
45 : G4VFastSimulationModel(aModelName, aEnvelope)
46 , fInference(new Par04InferenceSetup)
47 , fHitMaker(new G4FastSimHitMaker)
48 , fParallelHitMaker(new G4FastSimHitMaker) {
49 fParallelHitMaker->SetNameOfWorldWithSD("parallelWorldFastSim");
50}
51
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53
54Par04MLFastSimModel::Par04MLFastSimModel(G4String aModelName)
55 : G4VFastSimulationModel(aModelName)
56 , fInference(new Par04InferenceSetup)
57 , fHitMaker(new G4FastSimHitMaker)
58 , fParallelHitMaker(new G4FastSimHitMaker) {
59 fParallelHitMaker->SetNameOfWorldWithSD("parallelWorldFastSim");
60}
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63
64Par04MLFastSimModel::~Par04MLFastSimModel() {}
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67
68G4bool Par04MLFastSimModel::IsApplicable(const G4ParticleDefinition& aParticleType)
69{
70 return &aParticleType == G4Electron::ElectronDefinition() ||
71 &aParticleType == G4Positron::PositronDefinition() ||
72 &aParticleType == G4Gamma::GammaDefinition();
73}
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76
77G4bool Par04MLFastSimModel::ModelTrigger(const G4FastTrack& aFastTrack)
78{
79 return fInference->IfTrigger(aFastTrack.GetPrimaryTrack()->GetKineticEnergy());
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83
84void Par04MLFastSimModel::DoIt(const G4FastTrack& aFastTrack, G4FastStep& aFastStep)
85{
86 // remove particle from further processing by G4
87 aFastStep.KillPrimaryTrack();
88 aFastStep.SetPrimaryTrackPathLength(0.0);
89 G4double energy = aFastTrack.GetPrimaryTrack()->GetKineticEnergy();
90 aFastStep.SetTotalEnergyDeposited(energy);
91 G4ThreeVector position = aFastTrack.GetPrimaryTrack()->GetPosition();
92 G4ThreeVector direction = aFastTrack.GetPrimaryTrack()->GetMomentumDirection();
93
94 // calculate the incident angle
95 G4float angle = direction.theta();
96
97 // calculate how to deposit energy within the detector
98 // get it from inference model
99 fInference->GetEnergies(fEnergies, energy, angle);
100 fInference->GetPositions(fPositions, position, direction);
101
102 // deposit energy in the detector using calculated values of energy deposits
103 // and positions
104 for(size_t iHit = 0; iHit < fPositions.size(); iHit++)
105 {
106 if (fEnergies[iHit] > 0.0005) {
107 // Place hit in the physical readout of the detector
108 fParallelHitMaker->make(G4FastHit(fPositions[iHit], fEnergies[iHit]), aFastTrack);
109 // Place hit in the virtual grid, for validation purposes
110 fHitMaker->make(G4FastHit(fPositions[iHit], fEnergies[iHit]), aFastTrack);
111 }
112 }
113}
114#endif

Applications | User Support | Publications | Collaboration