source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/body/src/G4QCHIPSWorld.cc @ 836

Last change on this file since 836 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 5.6 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//
27// $Id: G4QCHIPSWorld.cc,v 1.32 2006/06/29 20:06:47 gunter Exp $
28// GEANT4 tag $Name:  $
29//
30//      ---------------- G4QCHIPSWorld ----------------
31//             by Mikhail Kossov, Sept 1999.
32//      class for the CHIPS World definition in CHIPS Model
33// -------------------------------------------------------------------
34
35//#define debug
36//#define pdebug
37
38#include "G4QCHIPSWorld.hh"
39
40// Initialization of the CHIPSWorld Pointer
41// G4QCHIPSWorld* G4QCHIPSWorld::aWorld =0;
42// G4QParticleVector G4QCHIPSWorld::qWorld;
43
44// Constructor
45G4QCHIPSWorld::G4QCHIPSWorld()
46{
47}
48
49G4QCHIPSWorld::~G4QCHIPSWorld()      // The CHIPS World is destructed only in the EndOfJob
50{
51  //G4int nP=GetQWorld().size();
52  //G4cout<<"G4QCHIPSWorld::Destructor: Before nP="<<nP<<","<<GetQWorld()[0]<<G4endl; //TMP
53  //if(nP) std::for_each(GetQWorld().begin(),GetQWorld().end(),DeleteQParticle());
54  //G4cout<<"G4QCHIPSWorld::Destructor: After"<<G4endl; // TMP
55  //GetQWorld().clear();
56}
57
58// Standard output for CHIPS World
59std::ostream& operator<<(std::ostream& lhs, G4QCHIPSWorld& rhs)
60//       ============================================
61{
62  // @@ Later make a list of activated particles and clusters
63  lhs << "[ Currently a#of particles in the CHIPS World = " << rhs.GetQPEntries() << "]";
64  return lhs;
65}
66
67// Returns Pointer to the CHIPS World
68G4QCHIPSWorld* G4QCHIPSWorld::Get()
69//             ====================
70{
71  static G4QCHIPSWorld theWorld;                // *** Static body of the CHIPS World ***
72// Returns Pointer to the CHIPS World
73//  if(!aWorld) aWorld=&theWorld;                 // Init the static pointer to CHIPS World
74//  return aWorld;
75  return &theWorld;
76}
77
78G4QParticleVector & G4QCHIPSWorld::GetQWorld()
79{
80  static G4QParticleVector theWorld;
81  return theWorld;
82}
83
84// Return pointer to Particles of the CHIPS World
85G4QParticleVector* G4QCHIPSWorld::GetParticles(G4int nOfParts)
86//                 ===========================================
87{
88  //static const G4int mnofParts = 486;           // max number of particles (up to A=80)
89  static const G4int mnofParts = 494;           // max number of particles (up to A=80) IN
90  static const G4bool cf = true;                // verbose=true G4QPDG construction flag
91#ifdef debug
92  G4cout<<"G4QCHIPSWorld::GetParticles: n="<<nOfParts<<" particles"<<G4endl;
93#endif
94  //if(GetQWorld().size())G4cout<<"G4QCHIPSWorld::GetPts:***Pt#0="<<GetQWorld()[0]<<G4endl;
95  if(nOfParts>0)
96  {
97#ifdef debug
98    G4cout<<"G4QCHIPSWorld::GetParticles: Creating CHIPS World of nP="<<nOfParts<<G4endl;
99#endif
100    G4int curNP=GetQWorld().size();
101    //G4cout<<"G4QCHIPSWorld::GetParticles: Creating CHIPS World of curNP="<<curNP<<G4endl;
102    if(nOfParts>curNP)                         // Initialization for increasing CHIPS World
103    {
104      if (nOfParts>mnofParts)
105      {
106        G4cerr<<"***G4QCHIPSWorld::GetPartics:nOfParts="<<nOfParts<<">"<<mnofParts<<G4endl;
107        nOfParts=mnofParts;
108      }
109      if (nOfParts<10) nOfParts=10;            // Minimal number of particles for Vacuum(?)
110#ifdef debug
111      G4cout<<"G4QCHIPSWorld::GetParticles: n="<<nOfParts<<",c="<<curNP<<G4endl;
112#endif
113      for (G4int i=curNP; i<nOfParts; i++) 
114      {
115#ifdef debug
116                G4cout<<"G4QCHIPSWorld::GetParticles: Create particle QCode="<<i<<G4endl;
117#endif
118        G4QParticle* curPart = new G4QParticle(cf,i); // Created with QCode=i
119#ifdef debug
120                G4cout<<"G4QCHIPSWorld::GetParticles: Particle QCode="<<i<<" is created."<<G4endl;
121#endif
122        //curPart->InitQParticle(i);             //
123        //if(!i) G4cout<<"G4QCHIPSWorld::GetParticles:Pt#0="<<curPart<<G4endl;
124        GetQWorld().push_back(curPart);           // Filled (forever but only once)
125#ifdef debug
126        G4cout<<"G4QCHIPSWorld::GetParticles: Pt#"<<i<<"("<<nOfParts<<") is done"<<G4endl;
127#endif
128      }
129    }
130    //else init--;//Recover theReInitializationCounter, if nothingWasAdded to theCHIPSWorld
131  }
132#ifdef debug
133  G4cout<<"G4QCHIPSWorld::GetParticles: TotalPt#"<<GetQWorld().size()<<G4endl;
134#endif
135  return &GetQWorld();
136}
137
138
139
140
141
Note: See TracBrowser for help on using the repository browser.