Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
B3DetectorConstruction Class Reference

#include <B3DetectorConstruction.hh>

Inheritance diagram for B3DetectorConstruction:
G4VUserDetectorConstruction

Public Member Functions

 B3DetectorConstruction ()
 
virtual ~B3DetectorConstruction ()
 
virtual G4VPhysicalVolumeConstruct ()
 
virtual void ConstructSDandField ()
 
- Public Member Functions inherited from G4VUserDetectorConstruction
 G4VUserDetectorConstruction ()
 
virtual ~G4VUserDetectorConstruction ()
 
virtual void CloneSD ()
 
virtual void CloneF ()
 
void RegisterParallelWorld (G4VUserParallelWorld *)
 
G4int ConstructParallelGeometries ()
 
void ConstructParallelSD ()
 
G4int GetNumberOfParallelWorld () const
 
G4VUserParallelWorldGetParallelWorld (G4int i) const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VUserDetectorConstruction
void SetSensitiveDetector (const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
 
void SetSensitiveDetector (G4LogicalVolume *logVol, G4VSensitiveDetector *aSD)
 

Detailed Description

Detector construction class to define materials and geometry.

Crystals are positioned in Ring, with an appropriate rotation matrix. Several copies of Ring are placed in the full detector.

Definition at line 45 of file B3DetectorConstruction.hh.

Constructor & Destructor Documentation

B3DetectorConstruction::B3DetectorConstruction ( )

Definition at line 51 of file B3DetectorConstruction.cc.

53  fCheckOverlaps(true)
54 {
55  DefineMaterials();
56 }
B3DetectorConstruction::~B3DetectorConstruction ( )
virtual

Definition at line 60 of file B3DetectorConstruction.cc.

61 { }

Member Function Documentation

G4VPhysicalVolume * B3DetectorConstruction::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Definition at line 83 of file B3DetectorConstruction.cc.

References python.hepunit::cm, python.hepunit::deg, G4NistManager::FindOrBuildMaterial(), G4cout, G4endl, G4Material::GetMaterialTable(), G4NistManager::Instance(), G4VisAttributes::Invisible, python.hepunit::mm, CLHEP::HepRotation::rotateY(), CLHEP::HepRotation::rotateZ(), G4LogicalVolume::SetVisAttributes(), and python.hepunit::twopi.

84 {
85  // Gamma detector Parameters
86  //
87  G4double cryst_dX = 6*cm, cryst_dY = 6*cm, cryst_dZ = 3*cm;
88  G4int nb_cryst = 32;
89  G4int nb_rings = 9;
90  //
91  G4double dPhi = twopi/nb_cryst, half_dPhi = 0.5*dPhi;
92  G4double cosdPhi = std::cos(half_dPhi);
93  G4double tandPhi = std::tan(half_dPhi);
94  //
95  G4double ring_R1 = 0.5*cryst_dY/tandPhi;
96  G4double ring_R2 = (ring_R1+cryst_dZ)/cosdPhi;
97  //
98  G4double detector_dZ = nb_rings*cryst_dX;
99  //
101  G4Material* default_mat = nist->FindOrBuildMaterial("G4_AIR");
102  G4Material* cryst_mat = nist->FindOrBuildMaterial("Lu2SiO5");
103 
104  //
105  // World
106  //
107  G4double world_sizeXY = 2.4*ring_R2;
108  G4double world_sizeZ = 1.2*detector_dZ;
109 
110  G4Box* solidWorld =
111  new G4Box("World", //its name
112  0.5*world_sizeXY, 0.5*world_sizeXY, 0.5*world_sizeZ); //its size
113 
114  G4LogicalVolume* logicWorld =
115  new G4LogicalVolume(solidWorld, //its solid
116  default_mat, //its material
117  "World"); //its name
118 
119  G4VPhysicalVolume* physWorld =
120  new G4PVPlacement(0, //no rotation
121  G4ThreeVector(), //at (0,0,0)
122  logicWorld, //its logical volume
123  "World", //its name
124  0, //its mother volume
125  false, //no boolean operation
126  0, //copy number
127  fCheckOverlaps); // checking overlaps
128 
129  //
130  // ring
131  //
132  G4Tubs* solidRing =
133  new G4Tubs("Ring", ring_R1, ring_R2, 0.5*cryst_dX, 0., twopi);
134 
135  G4LogicalVolume* logicRing =
136  new G4LogicalVolume(solidRing, //its solid
137  default_mat, //its material
138  "Ring"); //its name
139 
140  //
141  // define crystal
142  //
143  G4double gap = 0.5*mm; //a gap for wrapping
144  G4double dX = cryst_dX - gap, dY = cryst_dY - gap;
145  G4Box* solidCryst = new G4Box("crystal", dX/2, dY/2, cryst_dZ/2);
146 
147  G4LogicalVolume* logicCryst =
148  new G4LogicalVolume(solidCryst, //its solid
149  cryst_mat, //its material
150  "CrystalLV"); //its name
151 
152  // place crystals within a ring
153  //
154  for (G4int icrys = 0; icrys < nb_cryst ; icrys++) {
155  G4double phi = icrys*dPhi;
157  rotm.rotateY(90*deg);
158  rotm.rotateZ(phi);
159  G4ThreeVector uz = G4ThreeVector(std::cos(phi), std::sin(phi),0.);
160  G4ThreeVector position = (ring_R1+0.5*cryst_dZ)*uz;
161  G4Transform3D transform = G4Transform3D(rotm,position);
162 
163  new G4PVPlacement(transform, //rotation,position
164  logicCryst, //its logical volume
165  "crystal", //its name
166  logicRing, //its mother volume
167  false, //no boolean operation
168  icrys, //copy number
169  fCheckOverlaps); // checking overlaps
170  }
171 
172  //
173  // full detector
174  //
175  G4Tubs* solidDetector =
176  new G4Tubs("Detector", ring_R1, ring_R2, 0.5*detector_dZ, 0., twopi);
177 
178  G4LogicalVolume* logicDetector =
179  new G4LogicalVolume(solidDetector, //its solid
180  default_mat, //its material
181  "Detector"); //its name
182 
183  //
184  // place rings within detector
185  //
186  G4double OG = -0.5*(detector_dZ + cryst_dX);
187  for (G4int iring = 0; iring < nb_rings ; iring++) {
188  OG += cryst_dX;
189  new G4PVPlacement(0, //no rotation
190  G4ThreeVector(0,0,OG), //position
191  logicRing, //its logical volume
192  "ring", //its name
193  logicDetector, //its mother volume
194  false, //no boolean operation
195  iring, //copy number
196  fCheckOverlaps); // checking overlaps
197  }
198 
199  //
200  // place detector in world
201  //
202  new G4PVPlacement(0, //no rotation
203  G4ThreeVector(), //at (0,0,0)
204  logicDetector, //its logical volume
205  "Detector", //its name
206  logicWorld, //its mother volume
207  false, //no boolean operation
208  0, //copy number
209  fCheckOverlaps); // checking overlaps
210 
211  //
212  // patient
213  //
214  G4double patient_radius = 8*cm;
215  G4double patient_dZ = 10*cm;
216  G4Material* patient_mat = nist->FindOrBuildMaterial("G4_BRAIN_ICRP");
217 
218  G4Tubs* solidPatient =
219  new G4Tubs("Patient", 0., patient_radius, 0.5*patient_dZ, 0., twopi);
220 
221  G4LogicalVolume* logicPatient =
222  new G4LogicalVolume(solidPatient, //its solid
223  patient_mat, //its material
224  "PatientLV"); //its name
225 
226  //
227  // place patient in world
228  //
229  new G4PVPlacement(0, //no rotation
230  G4ThreeVector(), //at (0,0,0)
231  logicPatient, //its logical volume
232  "Patient", //its name
233  logicWorld, //its mother volume
234  false, //no boolean operation
235  0, //copy number
236  fCheckOverlaps); // checking overlaps
237 
238  // Visualization attributes
239  //
242 
243  // Print materials
245 
246  //always return the physical World
247  //
248  return physWorld;
249 }
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepRotation G4RotationMatrix
Definition: G4Box.hh:63
Definition: G4Tubs.hh:84
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:564
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
G4GLOB_DLL std::ostream G4cout
HepGeom::Transform3D G4Transform3D
static const G4VisAttributes Invisible
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:92
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetVisAttributes(const G4VisAttributes *pVA)
void B3DetectorConstruction::ConstructSDandField ( )
virtual

Reimplemented from G4VUserDetectorConstruction.

Definition at line 253 of file B3DetectorConstruction.cc.

References G4SDManager::GetSDMpointer(), G4MultiFunctionalDetector::RegisterPrimitive(), G4VUserDetectorConstruction::SetSensitiveDetector(), and G4SDManager::SetVerboseLevel().

254 {
256 
257  // declare crystal as a MultiFunctionalDetector scorer
258  //
260  G4VPrimitiveScorer* primitiv1 = new G4PSEnergyDeposit("edep");
261  cryst->RegisterPrimitive(primitiv1);
262  SetSensitiveDetector("CrystalLV",cryst);
263 
264  // declare patient as a MultiFunctionalDetector scorer
265  //
266  G4MultiFunctionalDetector* patient = new G4MultiFunctionalDetector("patient");
267  G4VPrimitiveScorer* primitiv2 = new G4PSDoseDeposit("dose");
268  patient->RegisterPrimitive(primitiv2);
269  SetSensitiveDetector("PatientLV",patient);
270 }
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
void SetVerboseLevel(G4int vl)
Definition: G4SDManager.hh:90
void SetSensitiveDetector(const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40

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