source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/interface/include/G4ProtonAntiProtonAtRestChips.hh @ 1007

Last change on this file since 1007 was 962, checked in by garnier, 15 years ago

update processes

File size: 4.9 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#ifndef G4ProtonAntiProtonAtRestChips_h
27#define G4ProtonAntiProtonAtRestChips_h
28
29#include "globals.hh"
30#include "G4VRestProcess.hh"
31#include "G4ParticleTable.hh"
32#include "G4Quasmon.hh"
33#include "G4QHadronVector.hh"
34#include "G4ParticleChange.hh"
35#include "G4LorentzVector.hh"
36#include "G4DynamicParticle.hh"
37#include "G4IonTable.hh"
38#include "G4Neutron.hh"
39#include "G4StopElementSelector.hh"
40#include "G4ChiralInvariantPhaseSpace.hh"
41#include "G4HadronicProcessType.hh"
42
43class G4ProtonAntiProtonAtRestChips : public G4VRestProcess
44{
45  private:
46  // hide assignment operator as private
47      G4ProtonAntiProtonAtRestChips& operator=(const G4ProtonAntiProtonAtRestChips &right);
48      G4ProtonAntiProtonAtRestChips(const G4ProtonAntiProtonAtRestChips& );
49   
50  public:
51 
52     G4ProtonAntiProtonAtRestChips(const G4String& processName ="AntiProtonAnnihilationAtRest")
53      : G4VRestProcess (processName, fHadronic) 
54     {
55       SetProcessSubType(fHadronAtRest);
56     }
57 
58    ~G4ProtonAntiProtonAtRestChips() {}
59
60     G4bool IsApplicable(const G4ParticleDefinition& aParticle)
61     {
62       return ( &aParticle == G4AntiProton::AntiProtonDefinition() );
63     }
64
65  // null physics table
66     void BuildPhysicsTable(const G4ParticleDefinition&){}
67
68     G4double AtRestGetPhysicalInteractionLength(const G4Track&track,
69                                                 G4ForceCondition*condition);
70
71  // zero mean lifetime
72     G4double GetMeanLifeTime(const G4Track& aTrack,
73                              G4ForceCondition* condition) {return 0.0;}
74
75     G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&); 
76
77  private:
78    G4ChiralInvariantPhaseSpace theModel;
79    G4StopElementSelector theSelector; // Assume identical laws as for muons
80};
81
82inline
83G4VParticleChange * G4ProtonAntiProtonAtRestChips::
84AtRestDoIt(const G4Track& aTrack, const G4Step&aStep)
85{
86  // Create target
87  G4Element * theTarget = theSelector.GetElement(aTrack.GetMaterial());
88  G4Nucleus aTargetNucleus(theTarget->GetN() ,theTarget->GetZ());
89
90  // Check model validity - note this will be a sub-branch in the ordinary stopping @@@@@@
91  // in the long haul. @@@@@@
92  if(aTrack.GetDynamicParticle()->GetDefinition() != G4AntiProton::AntiProton())
93  {
94    throw G4HadronicException(__FILE__, __LINE__, "Calling G4ProtonAntiProtonAtRestChips with particle other than p-bar!!!");
95  }
96  if(aTargetNucleus.GetZ() != 1)
97  {
98    throw G4HadronicException(__FILE__, __LINE__, "Calling G4ProtonAntiProtonAtRestChips for target other than Hydrogen!!!");
99  }
100 
101  // Call chips
102  return theModel.ApplyYourself(aTrack, aTargetNucleus);
103}
104
105G4double G4ProtonAntiProtonAtRestChips::
106AtRestGetPhysicalInteractionLength(const G4Track&track,
107                                   G4ForceCondition*condition)
108{
109  ResetNumberOfInteractionLengthLeft();
110  *condition = NotForced;
111  currentInteractionLength = GetMeanLifeTime(track, condition);
112#ifdef CHIPSdebug
113  if ((currentInteractionLength <0.0) || (verboseLevel>2))
114  {
115    G4cout << "G4ProtonAntiProtonAtRestChips::AtRestGetPhysicalInteractionLength ";
116    G4cout << "[ " << GetProcessName() << "]" <<G4endl;
117    track.GetDynamicParticle()->DumpInfo();
118    G4cout << " in Material  " << track.GetMaterial()->GetName() <<G4endl;
119    G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl;
120  }
121#endif
122  return theNumberOfInteractionLengthLeft * currentInteractionLength;
123}
124#endif
Note: See TracBrowser for help on using the repository browser.