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

Last change on this file since 836 was 819, checked in by garnier, 16 years ago

import all except CVS

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