Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4STRead Class Reference

#include <G4STRead.hh>

Public Member Functions

G4LogicalVolumeRead (const G4String &, G4Material *mediumMaterial, G4Material *solidMaterial)
 

Detailed Description

Definition at line 49 of file G4STRead.hh.

Member Function Documentation

◆ Read()

G4LogicalVolume * G4STRead::Read ( const G4String name,
G4Material mediumMaterial,
G4Material solidMaterial 
)

Definition at line 275 of file G4STRead.cc.

278{
279 if(mediumMaterial == nullptr)
280 {
281 G4Exception("G4STRead::Read()", "InvalidSetup", FatalException,
282 "Pointer to medium material is not valid!");
283 }
284 if(solidMaterial == nullptr)
285 {
286 G4Exception("G4STRead::Read()", "InvalidSetup", FatalException,
287 "Pointer to solid material is not valid!");
288 }
289
290 solid_material = solidMaterial;
291
292 world_box = new G4Box("TessellatedWorldBox", kInfinity, kInfinity, kInfinity);
293 // We don't know the extent of the world yet!
294 world_volume = new G4LogicalVolume(world_box, mediumMaterial,
295 "TessellatedWorldLV", 0, 0, 0);
296 world_extent = G4ThreeVector(0, 0, 0);
297
298 ReadGeom(name + ".geom");
299 ReadTree(name + ".tree");
300
301 // Now setting the world extent ...
302 //
303 if(world_box->GetXHalfLength() > world_extent.x())
304 {
305 world_box->SetXHalfLength(world_extent.x());
306 }
307 if(world_box->GetYHalfLength() > world_extent.y())
308 {
309 world_box->SetYHalfLength(world_extent.y());
310 }
311 if(world_box->GetZHalfLength() > world_extent.z())
312 {
313 world_box->SetZHalfLength(world_extent.z());
314 }
315
316 return world_volume;
317}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
CLHEP::Hep3Vector G4ThreeVector
double z() const
double x() const
double y() const
Definition: G4Box.hh:56
G4double GetYHalfLength() const
G4double GetZHalfLength() const
void SetZHalfLength(G4double dz)
Definition: G4Box.cc:172
G4double GetXHalfLength() const
void SetYHalfLength(G4double dy)
Definition: G4Box.cc:149
void SetXHalfLength(G4double dx)
Definition: G4Box.cc:125

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