Geant4  10.02.p03
TSDetectorConstruction Class Reference

#include <TSDetectorConstruction.hh>

Inheritance diagram for TSDetectorConstruction:
Collaboration diagram for TSDetectorConstruction:

Public Types

typedef std::map< G4String, G4Material * > MaterialCollection_t
 
typedef std::set< G4LogicalVolume * > ScoringVolumes_t
 

Public Member Functions

 TSDetectorConstruction ()
 
virtual ~TSDetectorConstruction ()
 
G4VPhysicalVolumeConstruct ()
 
const G4ThreeVectorGetWorldDimensions () const
 
const ScoringVolumes_tGetScoringVolumes () const
 
const G4StringGetMFDName () const
 
- 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
 

Static Public Member Functions

static TSDetectorConstructionInstance ()
 

Protected Member Functions

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

Private Attributes

G4VPhysicalVolumefWorldPhys
 
ScoringVolumes_t fScoringVolumes
 
G4String fWorldMaterialName
 
G4String fTargetMaterialName
 
G4String fCasingMaterialName
 
G4ThreeVector fWorldDim
 
G4ThreeVector fTargetDim
 
G4ThreeVector fTargetSections
 
G4String fMfdName
 

Static Private Attributes

static TSDetectorConstructionfgInstance = 0
 

Detailed Description

Definition at line 74 of file TSDetectorConstruction.hh.

Member Typedef Documentation

◆ MaterialCollection_t

◆ ScoringVolumes_t

Constructor & Destructor Documentation

◆ TSDetectorConstruction()

TSDetectorConstruction::TSDetectorConstruction ( )

Definition at line 91 of file TSDetectorConstruction.cc.

92 : fWorldPhys(0),
93  fWorldMaterialName("G4_Galactic"),
94  fTargetMaterialName("G4_B"),
95  fCasingMaterialName("G4_WATER"),
96  fWorldDim(G4ThreeVector(0.5*m, 0.5*m, 0.5*m)),
97  fTargetDim(G4ThreeVector(0.5*m, 0.5*m, 0.5*m)),
99  fMfdName("Target_MFD")
100 {
101  fgInstance = this;
102 }
CLHEP::Hep3Vector G4ThreeVector
static TSDetectorConstruction * fgInstance
G4VPhysicalVolume * fWorldPhys
static const double m
Definition: G4SIunits.hh:128

◆ ~TSDetectorConstruction()

TSDetectorConstruction::~TSDetectorConstruction ( )
virtual

Definition at line 106 of file TSDetectorConstruction.cc.

107 {
108  fgInstance = 0;
109 }
static TSDetectorConstruction * fgInstance

Member Function Documentation

◆ Construct()

G4VPhysicalVolume * TSDetectorConstruction::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Definition at line 113 of file TSDetectorConstruction.cc.

114 {
116 }
virtual MaterialCollection_t ConstructMaterials()
virtual G4VPhysicalVolume * ConstructWorld(const MaterialCollection_t &)
Here is the call graph for this function:

◆ ConstructMaterials()

TSDetectorConstruction::MaterialCollection_t TSDetectorConstruction::ConstructMaterials ( )
protectedvirtual

Definition at line 121 of file TSDetectorConstruction.cc.

122 {
123  MaterialCollection_t materials;
125 
126  materials["World"] = nist->FindOrBuildMaterial(fWorldMaterialName);
127  materials["Target"] = nist->FindOrBuildMaterial(fTargetMaterialName);
128  materials["Casing"] = nist->FindOrBuildMaterial(fCasingMaterialName);
129 
130  return materials;
131 }
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4NistManager * Instance()
std::map< G4String, G4Material * > MaterialCollection_t
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructSDandField()

void TSDetectorConstruction::ConstructSDandField ( )
protectedvirtual

Reimplemented from G4VUserDetectorConstruction.

Definition at line 276 of file TSDetectorConstruction.cc.

277 {
278  //------------------------------------------------//
279  // MultiFunctionalDetector //
280  //------------------------------------------------//
281  // Define MultiFunctionalDetector with name.
283 
284  G4VPrimitiveScorer* edep = new G4PSEnergyDeposit("EnergyDeposit");
285  MFDet->RegisterPrimitive(edep);
286  G4VPrimitiveScorer* nstep = new G4PSNofStep("NumberOfSteps");
287  MFDet->RegisterPrimitive(nstep);
288 
289  // add scoring volumes
290  for(auto ite : fScoringVolumes)
291  {
292  SetSensitiveDetector(ite, MFDet);
293  }
294 }
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
Double_t edep
void SetSensitiveDetector(const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructWorld()

G4VPhysicalVolume * TSDetectorConstruction::ConstructWorld ( const MaterialCollection_t materials)
protectedvirtual

Definition at line 136 of file TSDetectorConstruction.cc.

137 {
138  G4UserLimits* steplimit = new G4UserLimits(0.1*(fTargetDim.z()
139  /fTargetSections.z()));
140  G4bool check_overlap = true;
141 
142  G4Box* world_solid = new G4Box("World",
143  0.5*fWorldDim.x(),
144  0.5*fWorldDim.y(),
145  0.5*fWorldDim.z());
146  G4LogicalVolume* world_log = new G4LogicalVolume(world_solid,
147  materials.find("World")
148  ->second,
149  "World");
150  fWorldPhys = new G4PVPlacement(0,
151  G4ThreeVector(0.),
152  "World",
153  world_log,
154  0, false, 0, check_overlap);
155 
156  G4int nz = fTargetSections.z();
157  G4int ny = fTargetSections.y();
158  G4int nx = fTargetSections.x();
159 
160  // spacing between sections
164 
165  //G4cout << "World has dimensions : "
166  //<< G4BestUnit(fWorldDim, "Length") << G4endl;
167 
168  //------------------------------------------------------------------------//
169  // Set Visual Attributes
170  //------------------------------------------------------------------------//
171  G4VisAttributes* red = new G4VisAttributes(G4Color(1., 0., 0., 1.0));
172  G4VisAttributes* green = new G4VisAttributes(G4Color(0., 1., 0., 0.25));
173  G4VisAttributes* blue = new G4VisAttributes(G4Color(0., 0., 1., 0.1));
174  G4VisAttributes* white = new G4VisAttributes(G4Color(1., 1., 1., 1.));
175 
176  white->SetVisibility(true);
177  red->SetVisibility(true);
178  green->SetVisibility(true);
179  blue->SetVisibility(true);
180 
181  white->SetForceWireframe(true);
182  red->SetForceSolid(true);
183  green->SetForceSolid(true);
184  blue->SetForceSolid(true);
185 
186  world_log->SetVisAttributes(white);
187 
188  for(G4int k = 0; k < nz; ++k)
189  for(G4int j = 0; j < ny; ++j)
190  for(G4int i = 0; i < nx; ++i)
191  {
192  // displacement of section
193  G4double dx
194  = 0.5*sx + static_cast<G4double>(i)*sx - 0.5*fWorldDim.x();
195  G4double dy
196  = 0.5*sy + static_cast<G4double>(j)*sy - 0.5*fWorldDim.y();
197  G4double dz
198  = 0.5*sz + static_cast<G4double>(k)*sz - 0.5*fWorldDim.z();
199  G4ThreeVector td = G4ThreeVector(dx, dy, -dz);
200  // make unique name
201  std::stringstream ss_name;
202  ss_name << "Target_" << i << "_" << j << "_" << k;
203 
204  G4Box* target_solid = new G4Box(ss_name.str(),
205  0.5*sx,
206  0.5*sy,
207  0.5*sz);
208 
209 
210  G4Material* target_material = 0;
211  G4bool is_casing = true;
212 
213  if(j == 0 || j+1 == ny || i == 0 || i+1 == nx ||
214  (nz > 1 && (k == 0 || k+1 == nz)))
215  target_material = materials.find("Casing")->second;
216  else {
217  target_material = materials.find("Target")->second;
218  is_casing = false;
219  }
220 
221  G4LogicalVolume* target_log =
222  new G4LogicalVolume(target_solid,
223  target_material,
224  ss_name.str());
225 
226  target_log->SetUserLimits(steplimit);
227 
228  new G4PVPlacement(0,
229  td,
230  ss_name.str(),
231  target_log,
232  fWorldPhys,
233  true,
234  k*nx*ny + j*nx + i,
235  check_overlap);
236 
237  fScoringVolumes.insert(target_log);
238 
239  if(is_casing)
240  target_log->SetVisAttributes(blue);
241  else
242  {
243  // making a checkerboard for kicks...
244  G4bool even_z = (k%2 == 0) ? true : false;
245  G4bool even_y = (j%2 == 0) ? true : false;
246  G4bool even_x = (i%2 == 0) ? true : false;
247 
248  G4VisAttributes* theColor = nullptr;
249 
250  if((even_z))
251  {
252  if((even_y && even_x) || (!even_y && !even_x))
253  theColor = red;
254  else
255  theColor = green;
256  } else // ! even_z
257  {
258  if((!even_y && even_x) || (even_y && !even_x))
259  theColor = red;
260  else
261  theColor = green;
262  }
263 
264  target_log->SetVisAttributes(theColor);
265  }
266 
267  }
268 
269 
270 
271  return fWorldPhys;
272 }
void SetForceWireframe(G4bool)
Definition: test07.cc:36
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:64
void SetVisibility(G4bool)
void SetUserLimits(G4UserLimits *pULimits)
void SetForceSolid(G4bool)
Definition: test07.cc:36
int G4int
Definition: G4Types.hh:78
bool G4bool
Definition: G4Types.hh:79
G4Colour G4Color
Definition: G4Color.hh:42
G4VPhysicalVolume * fWorldPhys
double x() const
double y() const
double z() const
double G4double
Definition: G4Types.hh:76
void SetVisAttributes(const G4VisAttributes *pVA)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetMFDName()

const G4String& TSDetectorConstruction::GetMFDName ( ) const
inline

Definition at line 91 of file TSDetectorConstruction.hh.

91 { return fMfdName; }
Here is the call graph for this function:

◆ GetScoringVolumes()

const ScoringVolumes_t& TSDetectorConstruction::GetScoringVolumes ( ) const
inline

Definition at line 89 of file TSDetectorConstruction.hh.

90  { return fScoringVolumes; }

◆ GetWorldDimensions()

const G4ThreeVector& TSDetectorConstruction::GetWorldDimensions ( ) const
inline

Definition at line 88 of file TSDetectorConstruction.hh.

88 { return fWorldDim; }
Here is the caller graph for this function:

◆ Instance()

TSDetectorConstruction * TSDetectorConstruction::Instance ( void  )
static

Definition at line 84 of file TSDetectorConstruction.cc.

85 {
86  return fgInstance;
87 }
static TSDetectorConstruction * fgInstance
Here is the caller graph for this function:

Member Data Documentation

◆ fCasingMaterialName

G4String TSDetectorConstruction::fCasingMaterialName
private

Definition at line 104 of file TSDetectorConstruction.hh.

◆ fgInstance

TSDetectorConstruction * TSDetectorConstruction::fgInstance = 0
staticprivate

Definition at line 99 of file TSDetectorConstruction.hh.

◆ fMfdName

G4String TSDetectorConstruction::fMfdName
private

Definition at line 108 of file TSDetectorConstruction.hh.

◆ fScoringVolumes

ScoringVolumes_t TSDetectorConstruction::fScoringVolumes
private

Definition at line 101 of file TSDetectorConstruction.hh.

◆ fTargetDim

G4ThreeVector TSDetectorConstruction::fTargetDim
private

Definition at line 106 of file TSDetectorConstruction.hh.

◆ fTargetMaterialName

G4String TSDetectorConstruction::fTargetMaterialName
private

Definition at line 103 of file TSDetectorConstruction.hh.

◆ fTargetSections

G4ThreeVector TSDetectorConstruction::fTargetSections
private

Definition at line 107 of file TSDetectorConstruction.hh.

◆ fWorldDim

G4ThreeVector TSDetectorConstruction::fWorldDim
private

Definition at line 105 of file TSDetectorConstruction.hh.

◆ fWorldMaterialName

G4String TSDetectorConstruction::fWorldMaterialName
private

Definition at line 102 of file TSDetectorConstruction.hh.

◆ fWorldPhys

G4VPhysicalVolume* TSDetectorConstruction::fWorldPhys
private

Definition at line 100 of file TSDetectorConstruction.hh.


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