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