|
|
|||
File indexing completed on 2026-06-03 07:56:46
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 Par01PiModel.cc 0027 /// \brief Implementation of the Par01PiModel class 0028 0029 #include "Par01PiModel.hh" 0030 0031 #include "G4Gamma.hh" 0032 #include "G4PionMinus.hh" 0033 #include "G4PionPlus.hh" 0034 0035 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0036 0037 Par01PiModel::Par01PiModel(G4Region* anEnvelope) 0038 : G4VFastSimulationModel("Par01PiModel", anEnvelope) 0039 { 0040 ; 0041 } 0042 0043 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0044 0045 Par01PiModel::~Par01PiModel() 0046 { 0047 ; 0048 } 0049 0050 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0051 0052 G4bool Par01PiModel::IsApplicable(const G4ParticleDefinition& particleType) 0053 { 0054 return &particleType == G4PionMinus::PionMinusDefinition() 0055 || &particleType == G4PionPlus::PionPlusDefinition(); 0056 } 0057 0058 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0059 0060 G4bool Par01PiModel::ModelTrigger(const G4FastTrack& fastTrack) 0061 { 0062 //------------------------------------------------------------- 0063 // UserTrigger() method: method which has to decide if 0064 // the parameterisation has to be applied. 0065 // Here ModelTrigger() asks the user (ie you) a 0/1 answer. 0066 // 0067 // Note that quantities like the local/global position/direction etc.. 0068 // are available at this level via the fastTrack parameter (allowing 0069 // to check distance from boundaries, see below to allow the decision) 0070 //-------------------------------------------------------------- 0071 G4cout << "\nPar01PiModel::ModelTrigger() called:" << G4endl; 0072 G4cout << "--------------------------------" << G4endl; 0073 G4cout << "(particle is a " << fastTrack.GetPrimaryTrack()->GetDefinition()->GetParticleName() 0074 << " )\n" 0075 << G4endl; 0076 0077 // -- Examples of available informations: 0078 0079 // -- position: 0080 G4cout << " Track position: " << fastTrack.GetPrimaryTrack()->GetPosition() 0081 << "(global coord.)" << fastTrack.GetPrimaryTrackLocalPosition() << "(in envelope coord.)" 0082 << G4endl; 0083 0084 // -- direction: 0085 G4cout << " Track direction:" << fastTrack.GetPrimaryTrack()->GetMomentum().unit() 0086 << "(global coord.)" << fastTrack.GetPrimaryTrackLocalDirection() << "(in envelope coord.)" 0087 << G4endl; 0088 0089 return true; 0090 } 0091 0092 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0093 0094 void Par01PiModel::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep) 0095 //-------------------------------------------------------- 0096 // 0097 // User method to code the parameterisation properly said. 0098 // 0099 //-------------------------------------------------------- 0100 { 0101 //------------------------------------------------ 0102 // The primary track continues along its direction. 0103 // One secondary (a photon) is added: 0104 //------------------------------------------------ 0105 G4cout << " Pion `model' applied\n" << G4endl; 0106 0107 //------------------------------ 0108 // Primary: 0109 // idem as in "DefaultModel": 0110 // 0111 //------------------------------ 0112 G4ThreeVector position; 0113 G4double distance; 0114 distance = fastTrack.GetEnvelopeSolid()->DistanceToOut(fastTrack.GetPrimaryTrackLocalPosition(), 0115 fastTrack.GetPrimaryTrackLocalDirection()); 0116 position = 0117 fastTrack.GetPrimaryTrackLocalPosition() + distance * fastTrack.GetPrimaryTrackLocalDirection(); 0118 0119 // -- set final position: 0120 fastStep.ProposePrimaryTrackFinalPosition(position); 0121 0122 //--------------------------- 0123 // Secondary: 0124 // Adds one "secondary": 0125 //--------------------------- 0126 // -- First, user has to say how many secondaries will be created: 0127 fastStep.SetNumberOfSecondaryTracks(1); 0128 0129 //------------------------ 0130 // -- Build the secondary: 0131 //------------------------ 0132 // -- direction: 0133 G4ParticleMomentum direction(fastTrack.GetPrimaryTrackLocalDirection()); 0134 direction.setZ(direction.z() * 0.5); 0135 direction.setY(direction.y() + direction.z() * 0.1); 0136 direction = direction.unit(); // necessary ? 0137 0138 // -- dynamics (Note that many constructors exists for G4DynamicParticle 0139 // -- see prototype/particle+matter/particles/management/include/G4DynamicParticle.hh) 0140 G4DynamicParticle dynamique(G4Gamma::GammaDefinition(), direction, 0141 fastTrack.GetPrimaryTrack()->GetKineticEnergy() / 2.); 0142 // -- position: 0143 G4double Dist; 0144 Dist = fastTrack.GetEnvelopeSolid()->DistanceToOut(fastTrack.GetPrimaryTrackLocalPosition(), 0145 direction); 0146 G4ThreeVector posi; 0147 posi = fastTrack.GetPrimaryTrackLocalPosition() + Dist * direction; 0148 0149 //------------------------------------ 0150 //-- Creation of the secondary Track: 0151 //------------------------------------ 0152 fastStep.CreateSecondaryTrack(dynamique, posi, fastTrack.GetPrimaryTrack()->GetGlobalTime()); 0153 }
| [ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
|
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
|