Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Par02FastSimModelHCal Class Reference

#include <Par02FastSimModelHCal.hh>

Inheritance diagram for Par02FastSimModelHCal:
Collaboration diagram for Par02FastSimModelHCal:

Public Member Functions

 Par02FastSimModelHCal (G4String aModelName, G4Region *aEnvelope, Par02DetectorParametrisation::Parametrisation aParamType)
 
 Par02FastSimModelHCal (G4String aModelName, G4Region *aEnvelope)
 
 Par02FastSimModelHCal (G4String aModelName)
 
 ~Par02FastSimModelHCal ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &aParticle)
 
virtual G4bool ModelTrigger (const G4FastTrack &aFastTrack)
 
virtual void DoIt (const G4FastTrack &aFastTrack, G4FastStep &aFastStep)
 
- Public Member Functions inherited from G4VFastSimulationModel
 G4VFastSimulationModel (const G4String &aName)
 
 G4VFastSimulationModel (const G4String &aName, G4Envelope *, G4bool IsUnique=FALSE)
 
virtual ~G4VFastSimulationModel ()
 
virtual G4bool AtRestModelTrigger (const G4FastTrack &)
 
virtual void AtRestDoIt (const G4FastTrack &, G4FastStep &)
 
const G4String GetName () const
 
G4bool operator== (const G4VFastSimulationModel &) const
 

Detailed Description

Shortcut to the ordinary tracking for hadronic calorimeters.

Fast simulation model describes what should be done instead of a normal tracking. Instead of the ordinary tracking, a particle deposits its energy at the entrance to the hadronic calorimeter and its value is smeared (by Par02Smearer::SmearMomentum()). Based on G4 examples/extended/parametrisations/Par01/include/Par01EMShowerModel.hh .

Author
Anna Zaborowska

Definition at line 47 of file Par02FastSimModelHCal.hh.

Constructor & Destructor Documentation

Par02FastSimModelHCal::Par02FastSimModelHCal ( G4String  aModelName,
G4Region aEnvelope,
Par02DetectorParametrisation::Parametrisation  aParamType 
)

A constructor.

Parameters
aModelNameA name of the fast simulation model.
aEnvelopeA region where the model can take over the ordinary tracking.
aParamTypeA parametrisation type.

Definition at line 47 of file Par02FastSimModelHCal.cc.

48  :
49  G4VFastSimulationModel( aModelName, aEnvelope ), fCalculateParametrisation(),
50  fParametrisation( aType ) {}
G4VFastSimulationModel(const G4String &aName)
Par02FastSimModelHCal::Par02FastSimModelHCal ( G4String  aModelName,
G4Region aEnvelope 
)

A constructor.

Parameters
aModelNameA name of the fast simulation model.
aEnvelopeA region where the model can take over the ordinary tracking.

Definition at line 54 of file Par02FastSimModelHCal.cc.

55  :
56  G4VFastSimulationModel( aModelName, aEnvelope ), fCalculateParametrisation(),
57  fParametrisation( Par02DetectorParametrisation::eCMS ) {}
G4VFastSimulationModel(const G4String &aName)
Par02FastSimModelHCal::Par02FastSimModelHCal ( G4String  aModelName)

A constructor.

Parameters
aModelNameA name of the fast simulation model.

Definition at line 61 of file Par02FastSimModelHCal.cc.

61  :
62  G4VFastSimulationModel( aModelName ), fCalculateParametrisation(),
63  fParametrisation( Par02DetectorParametrisation::eCMS ) {}
G4VFastSimulationModel(const G4String &aName)
Par02FastSimModelHCal::~Par02FastSimModelHCal ( )

Definition at line 67 of file Par02FastSimModelHCal.cc.

67 {}

Member Function Documentation

void Par02FastSimModelHCal::DoIt ( const G4FastTrack aFastTrack,
G4FastStep aFastStep 
)
virtual

Smears the energy deposit and saves it, together with the position of the deposit, the hadronic calorimeter resolution and efficiency to the Par02PrimaryParticleInformation.

Parameters
aFastTrackA track.
aFastStepA step.

Implements G4VFastSimulationModel.

Definition at line 93 of file Par02FastSimModelHCal.cc.

94  {
95  //G4cout << " ________HCal model triggered _________" << G4endl;
96 
97  // Kill the parameterised particle at the entrance of the hadronic calorimeter
98  aFastStep.KillPrimaryTrack();
99  aFastStep.ProposePrimaryTrackPathLength( 0.0 );
100  G4double Edep = aFastTrack.GetPrimaryTrack()->GetKineticEnergy();
101 
102  // Consider only primary tracks (do nothing else for secondary hadrons)
103  G4ThreeVector Pos = aFastTrack.GetPrimaryTrack()->GetPosition();
104  if ( ! aFastTrack.GetPrimaryTrack()->GetParentID() ) {
107  if ( info->GetDoSmearing() ) {
108  // Smearing according to the hadronic calorimeter resolution
109  G4ThreeVector Porg = aFastTrack.GetPrimaryTrack()->GetMomentum();
110  G4double res = fCalculateParametrisation->
111  GetResolution( Par02DetectorParametrisation::eHCAL,
112  fParametrisation, Porg.mag() );
113  G4double eff = fCalculateParametrisation->
114  GetEfficiency( Par02DetectorParametrisation::eHCAL,
115  fParametrisation, Porg.mag() );
116  G4double Esm;
117  Esm = std::abs( Par02Smearer::Instance()->
118  SmearEnergy( aFastTrack.GetPrimaryTrack(), res ) );
119  Par02Output::Instance()->FillHistogram( 2, (Esm/MeV) / (Edep/MeV) );
120  // Setting the values of Pos, Esm, res and eff
121  Par02PrimaryParticleInformation* primaryInfo=
122  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  } else {
133  // No smearing: simply setting the value of Edep
134  ( (Par02PrimaryParticleInformation*) ( const_cast< G4PrimaryParticle* >
135  ( aFastTrack.GetPrimaryTrack()->GetDynamicParticle()->GetPrimaryParticle() )->
136  GetUserInformation() ) )->SetHCalEnergy( Edep );
137  // The (initial) energy of the particle is deposited in the step
138  // (which corresponds to the entrance of the hadronic calorimeter)
139  aFastStep.ProposeTotalEnergyDeposited( Edep );
140  }
141  }
142 }
const XML_Char XML_Encoding * info
Definition: expat.h:530
G4int GetParentID() const
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:208
static Par02Smearer * Instance()
Definition: Par02Smearer.cc:59
const G4DynamicParticle * GetDynamicParticle() const
ush Pos
Definition: deflate.h:89
const G4ThreeVector & GetPosition() const
void FillHistogram(G4int HNo, G4double value) const
Definition: Par02Output.cc:219
G4double GetKineticEnergy() const
void ProposePrimaryTrackPathLength(G4double)
static Par02Output * Instance()
Definition: Par02Output.cc:60
G4PrimaryParticle * GetPrimaryParticle() const
G4ThreeVector GetMomentum() const
static G4EventManager * GetEventManager()
void ProposeTotalEnergyDeposited(G4double anEnergyPart)
static constexpr double MeV
Definition: G4SIunits.hh:214
void KillPrimaryTrack()
Definition: G4FastStep.cc:88
double G4double
Definition: G4Types.hh:76
double mag() const
G4bool GetDoSmearing()
Gets the flag indicating if smearing should be done.
G4VUserEventInformation * GetUserInformation()

Here is the call graph for this function:

G4bool Par02FastSimModelHCal::IsApplicable ( const G4ParticleDefinition aParticle)
virtual

Checks if this model should be applied to this particle type.

Parameters
aParticleA particle definition (type).

Implements G4VFastSimulationModel.

Definition at line 71 of file Par02FastSimModelHCal.cc.

71  {
72  G4bool isOk = false;
73  // Applicable to all hadrons, i.e. any particle made of quarks
74  if ( aParticleType.GetQuarkContent(1) +
75  aParticleType.GetQuarkContent(2) +
76  aParticleType.GetQuarkContent(3) +
77  aParticleType.GetQuarkContent(4) +
78  aParticleType.GetQuarkContent(5) +
79  aParticleType.GetQuarkContent(6) != 0 ) {
80  isOk = true;
81  }
82  return isOk;
83 }
bool G4bool
Definition: G4Types.hh:79

Here is the call graph for this function:

G4bool Par02FastSimModelHCal::ModelTrigger ( const G4FastTrack aFastTrack)
virtual

Checks if the model should be applied, taking into account the kinematics of a track.

Parameters
aFastTrackA track.

Implements G4VFastSimulationModel.

Definition at line 87 of file Par02FastSimModelHCal.cc.

87  {
88  return true; // No kinematical restrictions to apply the parametrisation
89 }

The documentation for this class was generated from the following files: