source: trunk/examples/extended/biasing/B02/src/B02PSScoringDetectorConstruction.cc @ 1320

Last change on this file since 1320 was 1230, checked in by garnier, 14 years ago

update to geant4.9.3

File size: 8.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: B02PSScoringDetectorConstruction.cc,v 1.2 2007/06/22 13:38:55 ahoward Exp $
28// GEANT4 tag $Name: geant4-09-03-cand-01 $
29//
30
31#include "globals.hh"
32#include <sstream>
33
34#include "B02PSScoringDetectorConstruction.hh"
35
36#include "G4Material.hh"
37#include "G4Tubs.hh"
38#include "G4LogicalVolume.hh"
39#include "G4ThreeVector.hh"
40#include "G4PVPlacement.hh"
41
42// For Primitive Scorers
43#include "G4SDManager.hh"
44#include "G4MultiFunctionalDetector.hh"
45#include "G4SDParticleFilter.hh"
46#include "G4PSNofCollision.hh"
47#include "G4PSPopulation.hh"
48#include "G4PSTrackCounter.hh"
49#include "G4PSTrackLength.hh"
50
51
52
53B02PSScoringDetectorConstruction::B02PSScoringDetectorConstruction(G4String worldName)
54  :G4VUserParallelWorld(worldName),fLogicalVolumeVector()
55{
56  //  Construct();
57}
58
59B02PSScoringDetectorConstruction::~B02PSScoringDetectorConstruction()
60{
61  fLogicalVolumeVector.clear();
62}
63
64void B02PSScoringDetectorConstruction::Construct()
65{ 
66
67  G4cout << " constructing parallel world " << G4endl;
68
69  //GetWorld methods create a clone of the mass world to the parallel world (!)
70  // via the transportation manager
71  ghostWorld = GetWorld();
72  G4LogicalVolume* worldLogical = ghostWorld->GetLogicalVolume();
73  fLogicalVolumeVector.push_back(worldLogical);
74
75  G4String name("none");
76  G4double density(universe_mean_density), temperature(0), pressure(0);
77
78  name = "Galactic";
79  density     = universe_mean_density;            //from PhysicalConstants.h
80  pressure    = 3.e-18*pascal;
81  temperature = 2.73*kelvin;
82  G4cout << density << " " << kStateGas << G4endl;
83  G4Material *Galactic = 
84    new G4Material(name, 1., 1.01*g/mole, density,
85                   kStateGas,temperature,pressure);
86
87
88  //fPVolumeStore.AddPVolume(G4GeometryCell(*ghostWorld, 0));
89
90
91
92
93  // creating 18 slobs of 10 cm thicknes
94
95  G4double innerRadiusShield = 0*cm;
96  G4double outerRadiusShield = 100*cm;
97  G4double hightShield       = 5*cm;
98  G4double startAngleShield  = 0*deg;
99  G4double spanningAngleShield    = 360*deg;
100
101  G4Tubs *aShield = new G4Tubs("aShield",
102                               innerRadiusShield,
103                               outerRadiusShield,
104                               hightShield,
105                               startAngleShield,
106                               spanningAngleShield);
107 
108  // logical parallel cells
109
110  G4LogicalVolume *aShield_log = 
111    new G4LogicalVolume(aShield, Galactic, "aShield_log");
112  fLogicalVolumeVector.push_back(aShield_log);
113
114  // physical parallel cells
115
116  G4int i = 1;
117  G4double startz = -85*cm; 
118  for (i=1; i<=18; ++i) {
119   
120    name = GetCellName(i);
121   
122    G4double pos_x = 0*cm;
123    G4double pos_y = 0*cm;
124    G4double pos_z = startz + (i-1) * (2*hightShield);
125    G4VPhysicalVolume *pvol = 
126      new G4PVPlacement(0, 
127                        G4ThreeVector(pos_x, pos_y, pos_z),
128                        aShield_log, 
129                        name, 
130                        worldLogical, 
131                        false, 
132                        i);
133    //G4GeometryCell cell(*pvol, i);
134    //fPVolumeStore.AddPVolume(cell);
135  }
136
137  // filling the rest of the world volumr behind the concrete with
138  // another slob which should get the same importance value as the
139  // last slob
140  innerRadiusShield = 0*cm;
141  outerRadiusShield = 100*cm;
142  hightShield       = 5*cm;   
143  startAngleShield  = 0*deg;
144  spanningAngleShield    = 360*deg;
145
146  G4Tubs *aRest = new G4Tubs("Rest",
147                             innerRadiusShield,
148                             outerRadiusShield,
149                             hightShield,
150                             startAngleShield,
151                             spanningAngleShield);
152 
153  G4LogicalVolume *aRest_log = 
154    new G4LogicalVolume(aRest, Galactic, "aRest_log");
155  name = GetCellName(19);
156  fLogicalVolumeVector.push_back(aRest_log);
157   
158  G4double pos_x = 0*cm;
159  G4double pos_y = 0*cm;
160  G4double pos_z = 95*cm;
161  G4VPhysicalVolume *pvol = 
162    new G4PVPlacement(0, 
163                      G4ThreeVector(pos_x, pos_y, pos_z),
164                      aRest_log, 
165                      name, 
166                      worldLogical, 
167                      false, 
168                      19); // i = 19
169  //G4GeometryCell cell(*pvol, i);
170  //fPVolumeStore.AddPVolume(cell);
171
172  SetSensitive();
173
174}
175
176G4String B02PSScoringDetectorConstruction::GetCellName(G4int i) {
177  std::ostringstream os;
178  os << "cell_";
179  if (i<10) {
180    os << "0";
181  }
182  os << i;
183  G4String name = os.str();
184  return name;
185}
186
187G4VPhysicalVolume *B02PSScoringDetectorConstruction::GetWorldVolume() {
188   return ghostWorld;
189}
190
191
192G4VPhysicalVolume &B02PSScoringDetectorConstruction::GetWorldVolumeAddress() const{
193  return *ghostWorld;
194}
195
196void B02PSScoringDetectorConstruction::SetSensitive(){
197
198  //  -------------------------------------------------
199  //   The collection names of defined Primitives are
200  //   0       ConcreteSD/Collisions
201  //   1       ConcreteSD/CollWeight
202  //   2       ConcreteSD/Population
203  //   3       ConcreteSD/TrackEnter
204  //   4       ConcreteSD/SL
205  //   5       ConcreteSD/SLW
206  //   6       ConcreteSD/SLWE
207  //   7       ConcreteSD/SLW_V
208  //   8       ConcreteSD/SLWE_V
209  //  -------------------------------------------------
210
211
212  //================================================
213  // Sensitive detectors : MultiFunctionalDetector
214  //================================================
215  //
216  //  Sensitive Detector Manager.
217  G4SDManager* SDman = G4SDManager::GetSDMpointer();
218  //
219  // Sensitive Detector Name
220  G4String concreteSDname = "ConcreteSD";
221
222  //------------------------
223  // MultiFunctionalDetector
224  //------------------------
225  //
226  // Define MultiFunctionalDetector with name.
227  G4MultiFunctionalDetector* MFDet = new G4MultiFunctionalDetector(concreteSDname);
228  SDman->AddNewDetector( MFDet );                 // Register SD to SDManager
229
230
231  G4String fltName,particleName;
232  G4SDParticleFilter* neutronFilter = 
233      new G4SDParticleFilter(fltName="neutronFilter", particleName="neutron");
234
235  MFDet->SetFilter(neutronFilter);
236
237
238  for (std::vector<G4LogicalVolume *>::iterator it =  fLogicalVolumeVector.begin();
239       it != fLogicalVolumeVector.end(); it++){
240      (*it)->SetSensitiveDetector(MFDet);
241  }
242
243  G4String psName;
244  G4PSNofCollision*   scorer0 = new G4PSNofCollision(psName="Collisions"); 
245  MFDet->RegisterPrimitive(scorer0);
246
247
248  G4PSNofCollision*   scorer1 = new G4PSNofCollision(psName="CollWeight"); 
249  scorer1->Weighted(true);
250  MFDet->RegisterPrimitive(scorer1);
251
252
253  G4PSPopulation*   scorer2 = new G4PSPopulation(psName="Population"); 
254  MFDet->RegisterPrimitive(scorer2);
255
256  G4PSTrackCounter* scorer3 = new G4PSTrackCounter(psName="TrackEnter",fCurrent_In); 
257  MFDet->RegisterPrimitive(scorer3);
258
259  G4PSTrackLength* scorer4 = new G4PSTrackLength(psName="SL"); 
260  MFDet->RegisterPrimitive(scorer4);
261
262  G4PSTrackLength* scorer5 = new G4PSTrackLength(psName="SLW"); 
263  scorer5->Weighted(true);
264  MFDet->RegisterPrimitive(scorer5);
265
266  G4PSTrackLength* scorer6 = new G4PSTrackLength(psName="SLWE"); 
267  scorer6->Weighted(true);
268  scorer6->MultiplyKineticEnergy(true);
269  MFDet->RegisterPrimitive(scorer6);
270
271  G4PSTrackLength* scorer7 = new G4PSTrackLength(psName="SLW_V"); 
272  scorer7->Weighted(true);
273  scorer7->DivideByVelocity(true);
274  MFDet->RegisterPrimitive(scorer7);
275
276  G4PSTrackLength* scorer8 = new G4PSTrackLength(psName="SLWE_V"); 
277  scorer8->Weighted(true);
278  scorer8->MultiplyKineticEnergy(true);
279  scorer8->DivideByVelocity(true);
280  MFDet->RegisterPrimitive(scorer8);
281
282}
Note: See TracBrowser for help on using the repository browser.