source: trunk/source/processes/electromagnetic/xrays/include/G4Cerenkov.hh @ 1358

Last change on this file since 1358 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 8.5 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// $Id: G4Cerenkov.hh,v 1.11 2009/07/29 23:45:02 gum Exp $
28// GEANT4 tag $Name: geant4-09-04-beta-01 $
29//
30//
31////////////////////////////////////////////////////////////////////////
32// Cerenkov Radiation Class Definition
33////////////////////////////////////////////////////////////////////////
34//
35// File:        G4Cerenkov.hh 
36// Description: Discrete Process - Generation of Cerenkov Photons
37// Version:     2.0
38// Created:     1996-02-21
39// Author:      Juliet Armstrong
40// Updated:     2007-09-30 change inheritance to G4VDiscreteProcess
41//              2005-07-28 add G4ProcessType to constructor
42//              1999-10-29 add method and class descriptors
43//              1997-04-09 by Peter Gumplinger
44//              > G4MaterialPropertiesTable; new physics/tracking scheme
45// mail:        gum@triumf.ca
46//
47////////////////////////////////////////////////////////////////////////
48
49#ifndef G4Cerenkov_h
50#define G4Cerenkov_h 1
51
52/////////////
53// Includes
54/////////////
55
56#include "globals.hh"
57#include "templates.hh"
58#include "Randomize.hh"
59#include "G4ThreeVector.hh"
60#include "G4ParticleMomentum.hh"
61#include "G4Step.hh"
62#include "G4VProcess.hh"
63#include "G4OpticalPhoton.hh"
64#include "G4DynamicParticle.hh"
65#include "G4Material.hh"
66#include "G4PhysicsTable.hh"
67#include "G4MaterialPropertiesTable.hh"
68#include "G4PhysicsOrderedFreeVector.hh"
69
70// Class Description:
71// Discrete Process -- Generation of Cerenkov Photons.
72// Class inherits publicly from G4VDiscreteProcess.
73// Class Description - End:
74
75/////////////////////
76// Class Definition
77/////////////////////
78
79class G4Cerenkov : public G4VProcess
80{
81
82private:
83
84        //////////////
85        // Operators
86        //////////////
87
88        // G4Cerenkov& operator=(const G4Cerenkov &right);
89
90public: // Without description
91
92        ////////////////////////////////
93        // Constructors and Destructor
94        ////////////////////////////////
95
96        G4Cerenkov(const G4String& processName = "Cerenkov", 
97                            G4ProcessType type = fElectromagnetic);
98
99        // G4Cerenkov(const G4Cerenkov &right);
100
101        ~G4Cerenkov(); 
102
103        ////////////
104        // Methods
105        ////////////
106
107public: // With description
108
109        G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
110        // Returns true -> 'is applicable', for all charged particles
111        // except short-lived particles.
112
113        G4double GetMeanFreePath(const G4Track& aTrack,
114                                 G4double ,
115                                 G4ForceCondition* );
116        // Returns the discrete step limit and sets the 'StronglyForced'
117        // condition for the DoIt to be invoked at every step.
118
119        G4double PostStepGetPhysicalInteractionLength(const G4Track& aTrack,
120                                                      G4double ,
121                                                      G4ForceCondition* );
122        // Returns the discrete step limit and sets the 'StronglyForced'
123        // condition for the DoIt to be invoked at every step.
124
125        G4VParticleChange* PostStepDoIt(const G4Track& aTrack, 
126                                        const G4Step&  aStep);
127        // This is the method implementing the Cerenkov process.
128
129        //  no operation in  AtRestDoIt and  AlongStepDoIt
130        virtual G4double AlongStepGetPhysicalInteractionLength(
131                               const G4Track&,
132                               G4double  ,
133                               G4double  ,
134                               G4double& ,
135                               G4GPILSelection*
136                              ) { return -1.0; };
137
138        virtual G4double AtRestGetPhysicalInteractionLength(
139                               const G4Track& ,
140                               G4ForceCondition*
141                              ) { return -1.0; };
142
143        //  no operation in  AtRestDoIt and  AlongStepDoIt
144        virtual G4VParticleChange* AtRestDoIt(
145                               const G4Track& ,
146                               const G4Step&
147                              ) {return 0;};
148
149        virtual G4VParticleChange* AlongStepDoIt(
150                               const G4Track& ,
151                               const G4Step&
152                              ) {return 0;};
153
154        void SetTrackSecondariesFirst(const G4bool state);
155        // If set, the primary particle tracking is interrupted and any
156        // produced Cerenkov photons are tracked next. When all have
157        // been tracked, the tracking of the primary resumes.
158       
159        void SetMaxBetaChangePerStep(const G4double d);
160        // Set the maximum allowed change in beta = v/c in % (perCent)
161        // per step.
162
163        void SetMaxNumPhotonsPerStep(const G4int NumPhotons);
164        // Set the maximum number of Cerenkov photons allowed to be
165        // generated during a tracking step. This is an average ONLY;
166        // the actual number will vary around this average. If invoked,
167        // the maximum photon stack will roughly be of the size set.
168        // If not called, the step is not limited by the number of
169        // photons generated.
170
171        G4PhysicsTable* GetPhysicsTable() const;
172        // Returns the address of the physics table.
173
174        void DumpPhysicsTable() const;
175        // Prints the physics table.
176
177private:
178
179        void BuildThePhysicsTable();
180
181        /////////////////////
182        // Helper Functions
183        /////////////////////
184
185        G4double GetAverageNumberOfPhotons(const G4double charge,
186                                const G4double beta,
187                                const G4Material *aMaterial,
188                                const G4MaterialPropertyVector* Rindex) const;
189
190        ///////////////////////
191        // Class Data Members
192        ///////////////////////
193
194protected:
195
196        G4PhysicsTable* thePhysicsTable;
197        //  A Physics Table can be either a cross-sections table or
198        //  an energy table (or can be used for other specific
199        //  purposes).
200
201private:
202
203        G4bool fTrackSecondariesFirst;
204        G4double fMaxBetaChange;
205        G4int  fMaxPhotons;
206};
207
208////////////////////
209// Inline methods
210////////////////////
211
212inline 
213G4bool G4Cerenkov::IsApplicable(const G4ParticleDefinition& aParticleType)
214{
215   if (aParticleType.GetParticleName() == "chargedgeantino") return false;
216   if (aParticleType.IsShortLived()) return false;
217
218   return (aParticleType.GetPDGCharge() != 0);
219}
220
221inline 
222void G4Cerenkov::SetTrackSecondariesFirst(const G4bool state) 
223{ 
224        fTrackSecondariesFirst = state;
225}
226
227inline
228void G4Cerenkov::SetMaxBetaChangePerStep(const G4double value)
229{
230        fMaxBetaChange = value*perCent;
231}
232
233inline
234void G4Cerenkov::SetMaxNumPhotonsPerStep(const G4int NumPhotons) 
235{ 
236        fMaxPhotons = NumPhotons;
237}
238
239inline
240void G4Cerenkov::DumpPhysicsTable() const
241{
242        G4int PhysicsTableSize = thePhysicsTable->entries();
243        G4PhysicsOrderedFreeVector *v;
244
245        for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
246        {
247                v = (G4PhysicsOrderedFreeVector*)(*thePhysicsTable)[i];
248                v->DumpValues();
249        }
250}
251
252inline
253G4PhysicsTable* G4Cerenkov::GetPhysicsTable() const
254{
255  return thePhysicsTable;
256}
257
258#endif /* G4Cerenkov_h */
Note: See TracBrowser for help on using the repository browser.