Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4QCHIPSWorld.cc
Go to the documentation of this file.
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//
27// $Id$
28//
29// ---------------- G4QCHIPSWorld ----------------
30// by Mikhail Kossov, Sept 1999.
31// class for the CHIPS World definition in CHIPS Model
32// -------------------------------------------------------------------
33// Short description: The CHIPS World is a world of elementary particles
34// and nuclear fragments. This class is a singletone, but without fixed
35// limits. E.g. the nuclear fragments as possible G4Candidates can be
36// initialized in the CHIPS World only up to Be8 od C12 or other bigger
37// nuclear fragment. If one need the heavy fragment production then the
38// the CHIPS World must be initialized up to these fragments (see the
39// CHIPS Manual), but the price in performans will be big, because in
40// each act of the fragmentation competition these numerous candidates
41// take place in the competition and the hadronization probability is
42// calculated each time for each of them, so the Be8 limit (Be8->alpha+
43// alpha decays very fast and contribute to the alpha-spectrum) is the
44// most optimal.
45// -------------------------------------------------------------------
46
47//#define debug
48//#define pdebug
49
50#include "G4QCHIPSWorld.hh"
51
52// Initialization of the CHIPSWorld Pointer
53// G4QCHIPSWorld* G4QCHIPSWorld::aWorld =0;
54// G4QParticleVector G4QCHIPSWorld::qWorld;
55
56// Constructor
58{
59}
60
61G4QCHIPSWorld::~G4QCHIPSWorld() // The CHIPS World is destructed only in the EndOfJob
62{
63 //G4int nP=GetQWorld().size();
64 //G4cout<<"G4QCHIPSWorld::Destructor: Before nP="<<nP<<","<<GetQWorld()[0]<<G4endl; //TMP
65 //if(nP) std::for_each(GetQWorld().begin(),GetQWorld().end(),DeleteQParticle());
66 //G4cout<<"G4QCHIPSWorld::Destructor: After"<<G4endl; // TMP
67 //GetQWorld().clear();
68}
69
70// Standard output for CHIPS World
71std::ostream& operator<<(std::ostream& lhs, G4QCHIPSWorld& rhs)
72{
73 // @@ Later make a list of activated particles and clusters
74 lhs << "[ Currently a#of particles in the CHIPS World = " << rhs.GetQPEntries() << "]";
75 return lhs;
76}
77
78// Returns Pointer to the CHIPS World
80{
81 static G4QCHIPSWorld theWorld; // *** Static body of the CHIPS World ***
82// Returns Pointer to the CHIPS World
83// if(!aWorld) aWorld=&theWorld; // Init the static pointer to CHIPS World
84// return aWorld;
85 return &theWorld;
86}
87
88G4QParticleVector & G4QCHIPSWorld::GetQWorld()
89{
90 static G4QParticleVector theWorld;
91 return theWorld;
92}
93
94// Return pointer to Particles of the CHIPS World
96{
97 //static const G4int mnofParts = 486; // max number of particles (up to A=80)
98 //static const G4int mnofParts = 494; // max number of particles (up to A=80)
99 static const G4int mnofParts = 512; // max#of particles,A<57,G4QParticle::InitDecayVector
100 static const G4bool cf = true; // verbose=true G4QPDG construction flag
101#ifdef debug
102 G4cout<<"G4QCHIPSWorld::GetParticles: n="<<nOfParts<<" particles"<<G4endl;
103#endif
104 //if(GetQWorld().size())G4cout<<"G4QCHIPSWorld::GetPts:***Pt#0="<<GetQWorld()[0]<<G4endl;
105 if(nOfParts>0)
106 {
107#ifdef debug
108 G4cout<<"G4QCHIPSWorld::GetParticles: Creating CHIPS World of nP="<<nOfParts<<G4endl;
109#endif
110 G4int curNP=GetQWorld().size();
111 //G4cout<<"G4QCHIPSWorld::GetParticles: Creating CHIPS World of curNP="<<curNP<<G4endl;
112 if(nOfParts>curNP) // Initialization for increasing CHIPS World
113 {
114 if (nOfParts>mnofParts)
115 {
116 G4cerr<<"***G4QCHIPSWorld::GetPartics:nOfParts="<<nOfParts<<">"<<mnofParts<<G4endl;
117 nOfParts=mnofParts;
118 }
119 if (nOfParts<10) nOfParts=10; // Minimal number of particles for Vacuum(?)
120#ifdef debug
121 G4cout<<"G4QCHIPSWorld::GetParticles: n="<<nOfParts<<",c="<<curNP<<G4endl;
122#endif
123 for (G4int i=curNP; i<nOfParts; i++)
124 {
125#ifdef debug
126 G4cout<<"G4QCHIPSWorld::GetParticles: Create particle QCode="<<i<<G4endl;
127#endif
128 G4QParticle* curPart = new G4QParticle(cf,i); // Created with QCode=i
129#ifdef debug
130 G4cout<<"G4QCHIPSWorld::GetParticles: Particle QCode="<<i<<" is created."<<G4endl;
131#endif
132 //curPart->InitQParticle(i); //
133 //if(!i) G4cout<<"G4QCHIPSWorld::GetParticles:Pt#0="<<curPart<<G4endl;
134 GetQWorld().push_back(curPart); // Filled (forever but only once)
135#ifdef debug
136 G4cout<<"G4QCHIPSWorld::GetParticles: Pt#"<<i<<"("<<nOfParts<<") is done"<<G4endl;
137#endif
138 }
139 }
140 //else init--;//Recover theReInitializationCounter, if nothingWasAdded to theCHIPSWorld
141 }
142#ifdef debug
143 G4cout<<"G4QCHIPSWorld::GetParticles: TotalPt#"<<GetQWorld().size()<<G4endl;
144#endif
145 return &GetQWorld();
146}
std::ostream & operator<<(std::ostream &lhs, G4QCHIPSWorld &rhs)
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
G4DLLIMPORT std::ostream G4cout
static G4QCHIPSWorld * Get()
G4int GetQPEntries() const
G4QParticleVector * GetParticles(G4int nOfParts=0)