Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Private Attributes | List of all members
B03ImportanceDetectorConstruction Class Reference

#include <Doxymodules_biasing.h>

Inheritance diagram for B03ImportanceDetectorConstruction:
G4VUserParallelWorld

Public Member Functions

 B03ImportanceDetectorConstruction (G4String worldName)
 
 ~B03ImportanceDetectorConstruction ()
 
const G4VPhysicalVolumeGetPhysicalVolumeByName (const G4String &name) const
 
G4VPhysicalVolumeGetWorldVolumeAddress () const
 
G4String ListPhysNamesAsG4String ()
 
G4String GetCellName (G4int i)
 
G4GeometryCell GetGeometryCell (G4int i)
 
G4VPhysicalVolumeGetWorldVolume ()
 
void SetSensitive ()
 
virtual void ConstructSD ()
 

Private Member Functions

virtual void Construct ()
 

Private Attributes

B03PVolumeStore fPVolumeStore
 
std::vector< G4LogicalVolume * > fLogicalVolumeVector
 
G4VPhysicalVolumefGhostWorld
 

Detailed Description

Definition at line 51 of file Doxymodules_biasing.h.

Constructor & Destructor Documentation

◆ B03ImportanceDetectorConstruction()

B03ImportanceDetectorConstruction::B03ImportanceDetectorConstruction ( G4String  worldName)

Definition at line 56 of file B03ImportanceDetectorConstruction.cc.

59{
60 // Construct();
61}
std::vector< G4LogicalVolume * > fLogicalVolumeVector

◆ ~B03ImportanceDetectorConstruction()

B03ImportanceDetectorConstruction::~B03ImportanceDetectorConstruction ( )

Definition at line 65 of file B03ImportanceDetectorConstruction.cc.

66{
68}

Member Function Documentation

◆ GetPhysicalVolumeByName()

const G4VPhysicalVolume & B03ImportanceDetectorConstruction::GetPhysicalVolumeByName ( const G4String name) const

Definition at line 185 of file B03ImportanceDetectorConstruction.cc.

186 {
187 return *fPVolumeStore.GetPVolume(name);
188}
const G4VPhysicalVolume * GetPVolume(const G4String &name) const

◆ GetWorldVolumeAddress()

G4VPhysicalVolume & B03ImportanceDetectorConstruction::GetWorldVolumeAddress ( ) const

Definition at line 228 of file B03ImportanceDetectorConstruction.cc.

229 {
230 return *fGhostWorld;
231}

◆ ListPhysNamesAsG4String()

G4String B03ImportanceDetectorConstruction::ListPhysNamesAsG4String ( )

Definition at line 192 of file B03ImportanceDetectorConstruction.cc.

192 {
194 return names;
195}
G4String GetPNames() const

◆ GetCellName()

G4String B03ImportanceDetectorConstruction::GetCellName ( G4int  i)

Definition at line 199 of file B03ImportanceDetectorConstruction.cc.

199 {
200 std::ostringstream os;
201 os << "cell_";
202 if (i<10) {
203 os << "0";
204 }
205 os << i;
206 G4String name = os.str();
207 return name;
208}

◆ GetGeometryCell()

G4GeometryCell B03ImportanceDetectorConstruction::GetGeometryCell ( G4int  i)

Definition at line 212 of file B03ImportanceDetectorConstruction.cc.

212 {
213 G4String name(GetCellName(i));
214 const G4VPhysicalVolume *p=0;
215 p = fPVolumeStore.GetPVolume(name);
216 if (p) {
217 return G4GeometryCell(*p,0);
218 }
219 else {
220 G4cout << "B03ImportanceDetectorConstruction::GetGeometryCell: " << G4endl
221 << " couldn't get G4GeometryCell" << G4endl;
222 return G4GeometryCell(*fGhostWorld,-2);
223 }
224}

◆ GetWorldVolume()

G4VPhysicalVolume * B03ImportanceDetectorConstruction::GetWorldVolume ( )

Definition at line 235 of file B03ImportanceDetectorConstruction.cc.

235 {
236 return fGhostWorld;
237}

◆ SetSensitive()

void B03ImportanceDetectorConstruction::SetSensitive ( )

Definition at line 241 of file B03ImportanceDetectorConstruction.cc.

241 {
242
243 // -------------------------------------------------
244 // The collection names of defined Primitives are
245 // 0 ConcreteSD/Collisions
246 // 1 ConcreteSD/CollWeight
247 // 2 ConcreteSD/Population
248 // 3 ConcreteSD/TrackEnter
249 // 4 ConcreteSD/SL
250 // 5 ConcreteSD/SLW
251 // 6 ConcreteSD/SLWE
252 // 7 ConcreteSD/SLW_V
253 // 8 ConcreteSD/SLWE_V
254 // -------------------------------------------------
255
256 //now moved to ConstructSD()
257
258}

◆ ConstructSD()

void B03ImportanceDetectorConstruction::ConstructSD ( )
virtual

Definition at line 261 of file B03ImportanceDetectorConstruction.cc.

262{
263
264 G4SDManager* SDman = G4SDManager::GetSDMpointer();
265 //
266 // Sensitive Detector Name
267 G4String concreteSDname = "ConcreteSD";
268
269 //------------------------
270 // MultiFunctionalDetector
271 //------------------------
272 //
273 // Define MultiFunctionalDetector with name.
275 new G4MultiFunctionalDetector(concreteSDname);
276 SDman->AddNewDetector( MFDet ); // Register SD to SDManager
277
278 G4String fltName,particleName;
279 G4SDParticleFilter* neutronFilter =
280 new G4SDParticleFilter(fltName="neutronFilter", particleName="neutron");
281
282 MFDet->SetFilter(neutronFilter);
283
284 for (std::vector<G4LogicalVolume *>::iterator it =
285 fLogicalVolumeVector.begin();
286 it != fLogicalVolumeVector.end(); it++){
287 // (*it)->SetSensitiveDetector(MFDet);
288 SetSensitiveDetector((*it)->GetName(), MFDet);
289 }
290
291 G4String psName;
292 G4PSNofCollision* scorer0 = new G4PSNofCollision(psName="Collisions");
293 MFDet->RegisterPrimitive(scorer0);
294
295 G4PSNofCollision* scorer1 = new G4PSNofCollision(psName="CollWeight");
296 scorer1->Weighted(true);
297 MFDet->RegisterPrimitive(scorer1);
298
299 G4PSPopulation* scorer2 = new G4PSPopulation(psName="Population");
300 MFDet->RegisterPrimitive(scorer2);
301
302 G4PSTrackCounter* scorer3 =
303 new G4PSTrackCounter(psName="TrackEnter",fCurrent_In);
304 MFDet->RegisterPrimitive(scorer3);
305
306 G4PSTrackLength* scorer4 = new G4PSTrackLength(psName="SL");
307 MFDet->RegisterPrimitive(scorer4);
308
309 G4PSTrackLength* scorer5 = new G4PSTrackLength(psName="SLW");
310 scorer5->Weighted(true);
311 MFDet->RegisterPrimitive(scorer5);
312
313 G4PSTrackLength* scorer6 = new G4PSTrackLength(psName="SLWE");
314 scorer6->Weighted(true);
315 scorer6->MultiplyKineticEnergy(true);
316 MFDet->RegisterPrimitive(scorer6);
317
318 G4PSTrackLength* scorer7 = new G4PSTrackLength(psName="SLW_V");
319 scorer7->Weighted(true);
320 scorer7->DivideByVelocity(true);
321 MFDet->RegisterPrimitive(scorer7);
322
323 G4PSTrackLength* scorer8 = new G4PSTrackLength(psName="SLWE_V");
324 scorer8->Weighted(true);
325 scorer8->MultiplyKineticEnergy(true);
326 scorer8->DivideByVelocity(true);
327 MFDet->RegisterPrimitive(scorer8);
328}

◆ Construct()

void B03ImportanceDetectorConstruction::Construct ( )
privatevirtual

Definition at line 72 of file B03ImportanceDetectorConstruction.cc.

73{
74 G4cout << " constructing parallel world " << G4endl;
75
76 G4Material* dummyMat = 0;
77
78 //GetWorld methods create a clone of the mass world to the parallel world (!)
79 // via the transportation manager
80 fGhostWorld = GetWorld();
81 G4cout << " B03ImportanceDetectorConstruction:: ghostWorldName = "
82 << fGhostWorld->GetName() << G4endl;
83 G4LogicalVolume* worldLogical = fGhostWorld->GetLogicalVolume();
84 fLogicalVolumeVector.push_back(worldLogical);
85
86 G4String name("none");
87 // fPVolumeStore.AddPVolume(G4GeometryCell(*pWorldVolume, 0));
89
90 // creating 18 slobs of 10 cm thicknes
91
92 G4double innerRadiusShield = 0*cm;
93 G4double outerRadiusShield = 100*cm;
94 G4double heightShield = 5*cm;
95 G4double startAngleShield = 0*deg;
96 G4double spanningAngleShield = 360*deg;
97
98 G4Tubs *aShield = new G4Tubs("aShield",
99 innerRadiusShield,
100 outerRadiusShield,
101 heightShield,
102 startAngleShield,
103 spanningAngleShield);
104
105 // logical parallel cells
106
107 G4LogicalVolume *aShield_log_imp =
108 new G4LogicalVolume(aShield, dummyMat, "aShield_log_imp");
109 fLogicalVolumeVector.push_back(aShield_log_imp);
110
111 // physical parallel cells
112
113 G4int i = 1;
114 G4double startz = -85*cm;
115 // for (i=1; i<=18; ++i) {
116 for (i=1; i<=18; i++) {
117
118 name = GetCellName(i);
119
120 G4double pos_x = 0*cm;
121 G4double pos_y = 0*cm;
122 G4double pos_z = startz + (i-1) * (2*heightShield);
123 G4VPhysicalVolume *pvol =
124 new G4PVPlacement(0,
125 G4ThreeVector(pos_x, pos_y, pos_z),
126 aShield_log_imp,
127 name,
128 worldLogical,
129 false,
130 i);
131 // 0);
132 G4GeometryCell cell(*pvol, i);
133 // G4GeometryCell cell(*pvol, 0);
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 = 110*cm; exceeds world volume!!!!
142 outerRadiusShield = 100*cm;
143 // heightShield = 10*cm;
144 heightShield = 5*cm;
145 startAngleShield = 0*deg;
146 spanningAngleShield = 360*deg;
147
148 G4Tubs *aRest = new G4Tubs("Rest",
149 innerRadiusShield,
150 outerRadiusShield,
151 heightShield,
152 startAngleShield,
153 spanningAngleShield);
154
155 G4LogicalVolume *aRest_log =
156 new G4LogicalVolume(aRest, dummyMat, "aRest_log");
157
158 fLogicalVolumeVector.push_back(aRest_log);
159
160 name = GetCellName(19);
161
162 G4double pos_x = 0*cm;
163 G4double pos_y = 0*cm;
164 // G4double pos_z = 100*cm;
165 G4double pos_z = 95*cm;
166 G4VPhysicalVolume *pvol =
167 new G4PVPlacement(0,
168 G4ThreeVector(pos_x, pos_y, pos_z),
169 aRest_log,
170 name,
171 worldLogical,
172 false,
173 19);
174 // 0);
175 G4GeometryCell cell(*pvol, 19);
176 // G4GeometryCell cell(*pvol, 0);
178
179 SetSensitive();
180
181}
void AddPVolume(const G4GeometryCell &cell)

Member Data Documentation

◆ fPVolumeStore

B03PVolumeStore B03ImportanceDetectorConstruction::fPVolumeStore
private

Definition at line 66 of file B03ImportanceDetectorConstruction.hh.

◆ fLogicalVolumeVector

std::vector< G4LogicalVolume * > B03ImportanceDetectorConstruction::fLogicalVolumeVector
private

Definition at line 69 of file B03ImportanceDetectorConstruction.hh.

◆ fGhostWorld

G4VPhysicalVolume* B03ImportanceDetectorConstruction::fGhostWorld
private

Definition at line 73 of file B03ImportanceDetectorConstruction.hh.


The documentation for this class was generated from the following files:

Applications | User Support | Publications | Collaboration