![]() |
|
|||
File indexing completed on 2025-02-23 09:22:32
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 /// \file Par01/src/Par01PiModel.cc 0027 /// \brief Implementation of the Par01PiModel class 0028 // 0029 // 0030 // 0031 #include "Par01PiModel.hh" 0032 0033 #include "G4Gamma.hh" 0034 #include "G4PionMinus.hh" 0035 #include "G4PionPlus.hh" 0036 0037 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0038 0039 Par01PiModel::Par01PiModel(G4Region* anEnvelope) 0040 : G4VFastSimulationModel("Par01PiModel", anEnvelope) 0041 { 0042 ; 0043 } 0044 0045 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0046 0047 Par01PiModel::~Par01PiModel() 0048 { 0049 ; 0050 } 0051 0052 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0053 0054 G4bool Par01PiModel::IsApplicable(const G4ParticleDefinition& particleType) 0055 { 0056 return &particleType == G4PionMinus::PionMinusDefinition() 0057 || &particleType == G4PionPlus::PionPlusDefinition(); 0058 } 0059 0060 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0061 0062 G4bool Par01PiModel::ModelTrigger(const G4FastTrack& fastTrack) 0063 { 0064 //------------------------------------------------------------- 0065 // UserTrigger() method: method which has to decide if 0066 // the parameterisation has to be applied. 0067 // Here ModelTrigger() asks the user (ie you) a 0/1 answer. 0068 // 0069 // Note that quantities like the local/global position/direction etc.. 0070 // are available at this level via the fastTrack parameter (allowing 0071 // to check distance from boundaries, see below to allow the decision) 0072 //-------------------------------------------------------------- 0073 G4cout << "\nPar01PiModel::ModelTrigger() called:" << G4endl; 0074 G4cout << "--------------------------------" << G4endl; 0075 G4cout << "(particle is a " << fastTrack.GetPrimaryTrack()->GetDefinition()->GetParticleName() 0076 << " )\n" 0077 << G4endl; 0078 0079 // -- Examples of available informations: 0080 0081 // -- position: 0082 G4cout << " Track position: " << fastTrack.GetPrimaryTrack()->GetPosition() 0083 << "(global coord.)" << fastTrack.GetPrimaryTrackLocalPosition() << "(in envelope coord.)" 0084 << G4endl; 0085 0086 // -- direction: 0087 G4cout << " Track direction:" << fastTrack.GetPrimaryTrack()->GetMomentum().unit() 0088 << "(global coord.)" << fastTrack.GetPrimaryTrackLocalDirection() << "(in envelope coord.)" 0089 << G4endl; 0090 0091 return true; 0092 } 0093 0094 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0095 0096 void Par01PiModel::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep) 0097 //-------------------------------------------------------- 0098 // 0099 // User method to code the parameterisation properly said. 0100 // 0101 //-------------------------------------------------------- 0102 { 0103 //------------------------------------------------ 0104 // The primary track continues along its direction. 0105 // One secondary (a photon) is added: 0106 //------------------------------------------------ 0107 G4cout << " Pion `model' applied\n" << G4endl; 0108 0109 //------------------------------ 0110 // Primary: 0111 // idem as in "DefaultModel": 0112 // 0113 //------------------------------ 0114 G4ThreeVector position; 0115 G4double distance; 0116 distance = fastTrack.GetEnvelopeSolid()->DistanceToOut(fastTrack.GetPrimaryTrackLocalPosition(), 0117 fastTrack.GetPrimaryTrackLocalDirection()); 0118 position = 0119 fastTrack.GetPrimaryTrackLocalPosition() + distance * fastTrack.GetPrimaryTrackLocalDirection(); 0120 0121 // -- set final position: 0122 fastStep.ProposePrimaryTrackFinalPosition(position); 0123 0124 //--------------------------- 0125 // Secondary: 0126 // Adds one "secondary": 0127 //--------------------------- 0128 // -- First, user has to say how many secondaries will be created: 0129 fastStep.SetNumberOfSecondaryTracks(1); 0130 0131 //------------------------ 0132 // -- Build the secondary: 0133 //------------------------ 0134 // -- direction: 0135 G4ParticleMomentum direction(fastTrack.GetPrimaryTrackLocalDirection()); 0136 direction.setZ(direction.z() * 0.5); 0137 direction.setY(direction.y() + direction.z() * 0.1); 0138 direction = direction.unit(); // necessary ? 0139 0140 // -- dynamics (Note that many constructors exists for G4DynamicParticle 0141 // -- see prototype/particle+matter/particles/management/include/G4DynamicParticle.hh) 0142 G4DynamicParticle dynamique(G4Gamma::GammaDefinition(), direction, 0143 fastTrack.GetPrimaryTrack()->GetKineticEnergy() / 2.); 0144 // -- position: 0145 G4double Dist; 0146 Dist = fastTrack.GetEnvelopeSolid()->DistanceToOut(fastTrack.GetPrimaryTrackLocalPosition(), 0147 direction); 0148 G4ThreeVector posi; 0149 posi = fastTrack.GetPrimaryTrackLocalPosition() + Dist * direction; 0150 0151 //------------------------------------ 0152 //-- Creation of the secondary Track: 0153 //------------------------------------ 0154 fastStep.CreateSecondaryTrack(dynamique, posi, fastTrack.GetPrimaryTrack()->GetGlobalTime()); 0155 }
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
![]() ![]() |