source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/interface/src/G4QDiscProcessMixer.cc @ 1196

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

maj sur la beta de geant 4.9.3

File size: 5.7 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// $Id: G4QDiscProcessMixer.cc,v 1.8 2009/04/22 12:26:13 mkossov Exp $
27// GEANT4 tag $Name: geant4-09-03-beta-cand-01 $
28//
29//      ---------------- G4QDiscProcessMixer class -------------------
30//                 by Mikhail Kossov, Aug 2007.
31// G4QDiscProcessMixer class of the CHIPS Simulation Branch in GEANT4
32// ------------------------------------------------------------------------
33// Short description: universal mixer of processes (NOT models as in GHAD!)
34// depending on the application energy region (defined by users).
35// ------------------------------------------------------------------------
36
37//#define debug
38
39#include "G4QDiscProcessMixer.hh"
40
41// Constructor
42G4QDiscProcessMixer::G4QDiscProcessMixer(const G4String& name,
43                                         const G4ParticleDefinition* particle,
44                                         G4ProcessType pType):
45  G4VDiscreteProcess(name, pType), DPParticle(particle)
46{
47#ifdef debug
48  G4cout<<"G4QDiscProcessMixer::Constructor is called processName="<<name<<G4endl;
49#endif
50  if (verboseLevel>0) G4cout<<GetProcessName()<<" process is created "<<G4endl;
51}
52
53// Destructor
54G4QDiscProcessMixer::~G4QDiscProcessMixer()
55{
56  // Now the responsibility of deleting is deligated to the user, who created them
57  //for_each(theDPVector.begin(), theDPVector.end(), DeleteDiscreteProcess());
58}
59
60void G4QDiscProcessMixer::AddDiscreteProcess(G4VDiscreteProcess* DP, G4double ME)
61{
62  static const G4double maxEn = 1.E8*megaelectronvolt; // Conditional INF
63  if(!theDPVector.size()) // The first process in the DiscreteProcessVector (MaxE=INF)
64  {
65    std::pair<G4VDiscreteProcess*, G4double>* QDiscProc =
66      new std::pair<G4VDiscreteProcess*, G4double>(DP, maxEn);
67    theDPVector.push_back( QDiscProc );
68  }
69  else
70  {
71    if(ME<theDPVector[theDPVector.size()-1]->second)
72    {
73      std::pair<G4VDiscreteProcess*, G4double>* QDiscProc =
74        new std::pair<G4VDiscreteProcess*, G4double>(DP, ME);
75      theDPVector.push_back( QDiscProc );
76    }
77    else // Wrong Max Energy Order for the new process in the sequence of processes
78    {
79      G4cerr<<"G4QDiscProcessMixer::AddDiscreteProcess:LastMaxE("<<theDPVector.size()-1
80            <<")="<<theDPVector[theDPVector.size()-1]->second<<" <= MaxE="<<ME<<G4endl;
81      G4Exception("G4QDiscProcessMixer::AddDiscreteProcess: Wrong Max Energy Order");
82    }
83  }
84}
85
86G4bool G4QDiscProcessMixer::IsApplicable(const G4ParticleDefinition& particle) 
87{
88  if(particle == *DPParticle) return true;
89  return false;
90}
91
92G4double G4QDiscProcessMixer::PostStepGetPhysicalInteractionLength(const G4Track& Track,
93                                                                   G4double  PrevStSize,
94                                                                   G4ForceCondition* F)
95{
96  G4double kEn=Track.GetDynamicParticle()->GetKineticEnergy(); // Projectile kinetic energy
97  G4int maxDP=theDPVector.size();
98  if(maxDP) for(G4int i=maxDP-1; i>-1; i--) if(kEn < theDPVector[i]->second)
99    return theDPVector[i]->first->PostStepGetPhysicalInteractionLength(Track,PrevStSize,F);
100  return DBL_MAX;
101}
102
103// A fake class for compilation only
104G4double G4QDiscProcessMixer::GetMeanFreePath(const G4Track&, G4double, G4ForceCondition*)
105{
106  return DBL_MAX;
107}
108
109G4VParticleChange* G4QDiscProcessMixer::PostStepDoIt(const G4Track& Track,
110                                                     const G4Step& Step)
111{
112  G4double kEn=Track.GetDynamicParticle()->GetKineticEnergy(); // Projectile kinetic energy
113  G4int maxDP=theDPVector.size();
114  if(maxDP) for(G4int i=maxDP-1; i>-1; i--) if(kEn < theDPVector[i]->second)
115  {
116    //EnMomConservation= theDPVector[i]->first->GetEnegryMomentumConservation();
117    //nOfNeutrons      = theDPVector[i]->first->GetNumberOfNeutronsInTarget();
118    return theDPVector[i]->first->PostStepDoIt(Track, Step);
119  }
120  return G4VDiscreteProcess::PostStepDoIt(Track, Step);
121}
122
123//G4LorentzVector G4QDiscProcessMixer::GetEnegryMomentumConservation()
124//                                                              {return EnMomConservation;}
125
126//G4int G4QDiscProcessMixer::GetNumberOfNeutronsInTarget()
127//                                                                    {return nOfNeutrons;}
Note: See TracBrowser for help on using the repository browser.