source: trunk/source/geometry/magneticfield/test/NTST/src/NTSTLooperDeath.cc @ 1202

Last change on this file since 1202 was 1199, checked in by garnier, 15 years ago

nvx fichiers dans CVS

File size: 5.4 KB
Line 
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// -- Bogus -- BaBar Object-Oriented Geant-based Unified Simulation
28//
29// NTSTLooperDeath
30//
31// Description:
32//   This is a simple GEANT4 process that destroys any particle below
33//   the specified total kinetic energy that reverse direction (loops 180
34//   degress) in the x/y plane.
35//
36//   Based on BgsChargedLowEnergyDeath
37//
38// Author List:
39//   David Williams
40//
41// Modification History:
42//
43//-----------------------------------------------------------------------------
44
45#include "NTSTLooperDeath.hh"
46
47#include "G4TransportationManager.hh"
48#include "G4FieldManager.hh"
49#include "G4MagneticField.hh"
50
51//
52// Constructor
53//
54NTSTLooperDeath::NTSTLooperDeath( G4double theMinMomentum,
55                                  const char* name,
56                                  G4ProcessType type )
57  : G4VProcess( name, type ), minMomentum( theMinMomentum )
58{;}
59
60NTSTLooperDeath::~NTSTLooperDeath(){;}
61//
62// PostStepGetPhysicalInteractionLength
63//
64G4double NTSTLooperDeath::PostStepGetPhysicalInteractionLength( 
65                                                const G4Track& track,
66                                                G4double  , // previousStepSize,
67                                                G4ForceCondition* condition ) 
68{
69  const G4DynamicParticle *particle = track.GetDynamicParticle();
70       
71  *condition = NotForced;
72
73  //
74  // We don't touch any particle above the cut momentum
75  //
76  if (particle->GetTotalMomentum() > minMomentum) return DBL_MAX;
77 
78  //
79  // Nor do we touch any particle with small transverse component
80  // to their momentum. Here the cutoff is somewhat arbitrary.
81  //
82  G4double vperp = track.GetMomentumDirection().perp();
83  if (vperp < 0.1) return DBL_MAX;
84 
85  //
86  // How far in azimuthal angle do we need to go before the
87  // particle loops back on itself?
88  //
89  G4ThreeVector initialDir(track.GetVertexMomentumDirection());
90  G4ThreeVector dx = track.GetPosition() - track.GetVertexPosition();
91  G4double dxPerp = dx.perp();
92  if (dxPerp < 1E-6) return DBL_MAX;
93 
94  G4double dot = ( initialDir.x()*dx.x() + initialDir.y()*dx.y() )/dxPerp;
95 
96  if (dot < 0) return 0;                // Already done so
97 
98  G4double phi = pi - std::acos(dot);
99 
100  if (phi < 0) return 0;
101
102  //
103  // What is the radius of curvature?
104  //
105  // Only use the z component of the field to calculate this.
106  //
107  G4FieldManager *fieldManager = G4TransportationManager::GetTransportationManager()->GetFieldManager();
108 
109  if (!fieldManager->DoesFieldExist()) return DBL_MAX;
110 
111  //
112  // Assume the field is a magnetic field (i.e. returns
113  // a vector of three doubles of unit Telsa).
114  //
115  // There is no good way I know to confirm this, so it is purely
116  // a matter of faith. *** BEWARE ***
117  //
118  G4MagneticField *field = (G4MagneticField *)fieldManager->GetDetectorField();
119  G4double b[3];
120  G4ThreeVector pos(track.GetPosition());
121  G4double posv[3] = { pos.x(), pos.y(), pos.z() };
122  field->GetFieldValue( posv, b );
123 
124  //
125  // No field? Forget it!
126  //
127  if (std::fabs(b[2]) < 0.00001) return DBL_MAX;
128 
129  //
130  // Calculate radius of curvature, the usual way.
131  // Note G4 default units: mm, Telsa, MeV
132  //
133  // G4 suggestion: the constant below should be added to
134  // the geant4 list.
135  //
136  G4double radius = std::fabs(track.GetMomentum().perp()/299.79251/b[2]);
137 
138  //
139  // Convert this to a distance
140  //
141  return radius*phi/vperp;
142}
143
144
145//
146// PostStepDoit
147//
148G4VParticleChange *
149NTSTLooperDeath::PostStepDoIt( const G4Track &track, const G4Step & ) // step )
150{
151  pParticleChange->Initialize(track);
152
153  //
154  // This is a tough one. What should happen when we kill off a looper?
155  // For now: just deposit remaining energy (including mass).
156  //
157  const G4DynamicParticle *particle = track.GetDynamicParticle();
158  G4double energyDeposited = particle->GetTotalEnergy();
159
160  pParticleChange->ProposeTrackStatus( fStopAndKill );
161  pParticleChange->SetNumberOfSecondaries( 0 );
162  pParticleChange->ProposeLocalEnergyDeposit( energyDeposited );
163  ClearNumberOfInteractionLengthLeft();
164
165  return pParticleChange;
166}
167
168
169
170
171
172
Note: See TracBrowser for help on using the repository browser.