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

Last change on this file since 1340 was 1196, checked in by garnier, 15 years ago

update CVS release candidate geant4.9.3.01

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