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

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

update processes

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