source: trunk/examples/extended/biasing/B02/src/B02ImportanceDetectorConstruction.cc@ 1036

Last change on this file since 1036 was 807, checked in by garnier, 17 years ago

update

File size: 9.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: B02ImportanceDetectorConstruction.cc,v 1.11 2007/06/22 13:38:55 ahoward Exp $
28// GEANT4 tag $Name: $
29//
30
31#include "globals.hh"
32#include <sstream>
33
34#include "B02ImportanceDetectorConstruction.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
53B02ImportanceDetectorConstruction::B02ImportanceDetectorConstruction(G4String worldName)
54:G4VUserParallelWorld(worldName),fLogicalVolumeVector()
55{
56 // Construct();
57}
58
59B02ImportanceDetectorConstruction::~B02ImportanceDetectorConstruction()
60{
61 fLogicalVolumeVector.clear();
62}
63
64void B02ImportanceDetectorConstruction::Construct()
65{
66 G4cout << " constructing parallel world " << G4endl;
67
68 //GetWorld methods create a clone of the mass world to the parallel world (!)
69 // via the transportation manager
70 ghostWorld = GetWorld();
71 G4LogicalVolume* worldLogical = ghostWorld->GetLogicalVolume();
72 fLogicalVolumeVector.push_back(worldLogical);
73
74 G4String name("none");
75 G4double density(universe_mean_density), temperature(0), pressure(0);
76
77 name = "Galactic";
78 density = universe_mean_density; //from PhysicalConstants.h
79 pressure = 3.e-18*pascal;
80 temperature = 2.73*kelvin;
81 G4cout << density << " " << kStateGas << G4endl;
82 G4Material *Galactic =
83 new G4Material(name, 1., 1.01*g/mole, density,
84 kStateGas,temperature,pressure);
85
86
87 // fPVolumeStore.AddPVolume(G4GeometryCell(*pWorldVolume, 0));
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 for (i=1; i<=18; i++) {
120
121 name = GetCellName(i);
122
123 G4double pos_x = 0*cm;
124 G4double pos_y = 0*cm;
125 G4double pos_z = startz + (i-1) * (2*hightShield);
126 G4VPhysicalVolume *pvol =
127 new G4PVPlacement(0,
128 G4ThreeVector(pos_x, pos_y, pos_z),
129 aShield_log,
130 name,
131 worldLogical,
132 false,
133 i);
134 // 0);
135 G4GeometryCell cell(*pvol, i);
136 // G4GeometryCell cell(*pvol, 0);
137 fPVolumeStore.AddPVolume(cell);
138 }
139
140 // filling the rest of the world volumr behind the concrete with
141 // another slob which should get the same importance value as the
142 // last slob
143 innerRadiusShield = 0*cm;
144 // outerRadiusShield = 110*cm; exceeds world volume!!!!
145 outerRadiusShield = 101*cm;
146 // hightShield = 10*cm;
147 hightShield = 5*cm;
148 startAngleShield = 0*deg;
149 spanningAngleShield = 360*deg;
150
151 G4Tubs *aRest = new G4Tubs("Rest",
152 innerRadiusShield,
153 outerRadiusShield,
154 hightShield,
155 startAngleShield,
156 spanningAngleShield);
157
158 G4LogicalVolume *aRest_log =
159 new G4LogicalVolume(aRest, Galactic, "aRest_log");
160
161 fLogicalVolumeVector.push_back(aRest_log);
162
163 name = GetCellName(19);
164
165 G4double pos_x = 0*cm;
166 G4double pos_y = 0*cm;
167 // G4double pos_z = 100*cm;
168 G4double pos_z = 95*cm;
169 G4VPhysicalVolume *pvol =
170 new G4PVPlacement(0,
171 G4ThreeVector(pos_x, pos_y, pos_z),
172 aRest_log,
173 name,
174 worldLogical,
175 false,
176 19);
177 // 0);
178 G4GeometryCell cell(*pvol, 19);
179 // G4GeometryCell cell(*pvol, 0);
180 fPVolumeStore.AddPVolume(cell);
181
182 SetSensitive();
183
184}
185
186const G4VPhysicalVolume &B02ImportanceDetectorConstruction::
187GetPhysicalVolumeByName(const G4String& name) const {
188 return *fPVolumeStore.GetPVolume(name);
189}
190
191
192G4String B02ImportanceDetectorConstruction::ListPhysNamesAsG4String(){
193 G4String names(fPVolumeStore.GetPNames());
194 return names;
195}
196
197
198G4String B02ImportanceDetectorConstruction::GetCellName(G4int i) {
199 std::ostringstream os;
200 os << "cell_";
201 if (i<10) {
202 os << "0";
203 }
204 os << i;
205 G4String name = os.str();
206 return name;
207}
208
209G4GeometryCell B02ImportanceDetectorConstruction::GetGeometryCell(G4int i){
210 G4String name(GetCellName(i));
211 const G4VPhysicalVolume *p=0;
212 p = fPVolumeStore.GetPVolume(name);
213 if (p) {
214 return G4GeometryCell(*p,0);
215 }
216 else {
217 G4cout << "B02ImportanceDetectorConstruction::GetGeometryCell: couldn't get G4GeometryCell" << G4endl;
218 return G4GeometryCell(*ghostWorld,-2);
219 }
220}
221
222
223G4VPhysicalVolume &B02ImportanceDetectorConstruction::GetWorldVolumeAddress() const{
224 return *ghostWorld;
225}
226
227G4VPhysicalVolume *B02ImportanceDetectorConstruction::GetWorldVolume() {
228 return ghostWorld;
229}
230
231
232void B02ImportanceDetectorConstruction::SetSensitive(){
233
234 // -------------------------------------------------
235 // The collection names of defined Primitives are
236 // 0 ConcreteSD/Collisions
237 // 1 ConcreteSD/CollWeight
238 // 2 ConcreteSD/Population
239 // 3 ConcreteSD/TrackEnter
240 // 4 ConcreteSD/SL
241 // 5 ConcreteSD/SLW
242 // 6 ConcreteSD/SLWE
243 // 7 ConcreteSD/SLW_V
244 // 8 ConcreteSD/SLWE_V
245 // -------------------------------------------------
246
247
248 //================================================
249 // Sensitive detectors : MultiFunctionalDetector
250 //================================================
251 //
252 // Sensitive Detector Manager.
253
254 G4SDManager* SDman = G4SDManager::GetSDMpointer();
255 //
256 // Sensitive Detector Name
257 G4String concreteSDname = "ConcreteSD";
258
259 //------------------------
260 // MultiFunctionalDetector
261 //------------------------
262 //
263 // Define MultiFunctionalDetector with name.
264 G4MultiFunctionalDetector* MFDet = new G4MultiFunctionalDetector(concreteSDname);
265 SDman->AddNewDetector( MFDet ); // Register SD to SDManager
266
267
268 G4String fltName,particleName;
269 G4SDParticleFilter* neutronFilter =
270 new G4SDParticleFilter(fltName="neutronFilter", particleName="neutron");
271
272 MFDet->SetFilter(neutronFilter);
273
274
275 for (std::vector<G4LogicalVolume *>::iterator it = fLogicalVolumeVector.begin();
276 it != fLogicalVolumeVector.end(); it++){
277 (*it)->SetSensitiveDetector(MFDet);
278 }
279
280 G4String psName;
281 G4PSNofCollision* scorer0 = new G4PSNofCollision(psName="Collisions");
282 MFDet->RegisterPrimitive(scorer0);
283
284
285 G4PSNofCollision* scorer1 = new G4PSNofCollision(psName="CollWeight");
286 scorer1->Weighted(true);
287 MFDet->RegisterPrimitive(scorer1);
288
289
290 G4PSPopulation* scorer2 = new G4PSPopulation(psName="Population");
291 MFDet->RegisterPrimitive(scorer2);
292
293 G4PSTrackCounter* scorer3 = new G4PSTrackCounter(psName="TrackEnter",fCurrent_In);
294 MFDet->RegisterPrimitive(scorer3);
295
296 G4PSTrackLength* scorer4 = new G4PSTrackLength(psName="SL");
297 MFDet->RegisterPrimitive(scorer4);
298
299 G4PSTrackLength* scorer5 = new G4PSTrackLength(psName="SLW");
300 scorer5->Weighted(true);
301 MFDet->RegisterPrimitive(scorer5);
302
303 G4PSTrackLength* scorer6 = new G4PSTrackLength(psName="SLWE");
304 scorer6->Weighted(true);
305 scorer6->MultiplyKineticEnergy(true);
306 MFDet->RegisterPrimitive(scorer6);
307
308 G4PSTrackLength* scorer7 = new G4PSTrackLength(psName="SLW_V");
309 scorer7->Weighted(true);
310 scorer7->DivideByVelocity(true);
311 MFDet->RegisterPrimitive(scorer7);
312
313 G4PSTrackLength* scorer8 = new G4PSTrackLength(psName="SLWE_V");
314 scorer8->Weighted(true);
315 scorer8->MultiplyKineticEnergy(true);
316 scorer8->DivideByVelocity(true);
317 MFDet->RegisterPrimitive(scorer8);
318
319}
Note: See TracBrowser for help on using the repository browser.