#include <VisExamples/VisG4ExN02/src/VisG4ExampleDetector.h>
Public Member Functions | |
virtual G4VPhysicalVolume * | Construct (void) |
double | fullLength (void) const |
VisG4ExampleDetector (void) | |
~VisG4ExampleDetector (void) | |
Private Attributes | |
G4LogicalVolume * | m_chamberLV |
G4VPhysicalVolume * | m_chamberPV |
G4Box * | m_chamberSol |
VisG4ExampleField * | m_field |
G4Material * | m_matPb |
G4Material * | m_matXenon |
G4LogicalVolume * | m_targetLV |
G4VPhysicalVolume * | m_targetPV |
G4Box * | m_targetSol |
G4LogicalVolume * | m_trackerLV |
G4VPhysicalVolume * | m_trackerPV |
G4Box * | m_trackerSol |
G4LogicalVolume * | m_worldLV |
G4VPhysicalVolume * | m_worldPV |
G4Box * | m_worldSol |
Static Private Attributes | |
static const double | CHAMBER_SIZE = 20 * cm |
static const double | CHAMBER_SPACING = 80 * cm |
static const int | N_CHAMBERS = 5 |
static const double | TARGET_SIZE = 5 * cm |
static const double | TRACKER_SIZE = (N_CHAMBERS+1) * CHAMBER_SPACING |
static const double | WORLD_SIZE = 1.2 * (TARGET_SIZE + TRACKER_SIZE) |
(FIXME: Use Andrea's BarrelSaggin detector.)
Definition at line 25 of file VisG4ExampleDetector.h.
VisG4ExampleDetector::VisG4ExampleDetector | ( | void | ) |
Definition at line 35 of file VisG4ExampleDetector.cc.
00036 : m_worldSol (0), m_worldLV (0), m_worldPV (0), 00037 m_targetSol (0), m_targetLV (0), m_targetPV (0), 00038 m_trackerSol (0), m_trackerLV (0), m_trackerPV (0), 00039 m_chamberSol (0), m_chamberLV (0), m_chamberPV (0), 00040 m_matPb (0), m_matXenon (0), 00041 m_field (new VisG4ExampleField (G4ThreeVector (0, 0, 1.*tesla))) 00042 {}
VisG4ExampleDetector::~VisG4ExampleDetector | ( | void | ) |
Definition at line 44 of file VisG4ExampleDetector.cc.
References m_field.
00045 { delete m_field; }
G4VPhysicalVolume * VisG4ExampleDetector::Construct | ( | void | ) | [virtual] |
Definition at line 48 of file VisG4ExampleDetector.cc.
References CHAMBER_SIZE, CHAMBER_SPACING, g, m_chamberLV, m_chamberPV, m_chamberSol, m_matPb, m_matXenon, m_targetLV, m_targetPV, m_targetSol, m_trackerLV, m_trackerPV, m_trackerSol, m_worldLV, m_worldPV, m_worldSol, N_CHAMBERS, TARGET_SIZE, TRACKER_SIZE, WORLD_SIZE, and cmsRelvalreport::yellow().
00049 { 00050 // MATERIAL: Air, Pb, Xenon gas 00051 G4Material *air = new G4Material ("Air", 1.29*mg/cm3, 2); 00052 air->AddElement (new G4Element ("Nitrogen", "N", 7., 14.01*g/mole), .7); 00053 air->AddElement (new G4Element ("Oxygen", "O", 8., 16.00*g/mole), .3); 00054 00055 m_matPb = new G4Material ("Pb", 82., 207.19*g/mole, 11.35*g/cm3); 00056 m_matXenon = new G4Material ("XenonGas", 54., 131.29*g/mole, 5.458*mg/cm3, 00057 kStateGas, 293.15*kelvin, 1*atmosphere); 00058 00059 // SOLIDS, LVs and PVs 00060 // - World 00061 m_worldSol = new G4Box ("world", WORLD_SIZE/2, WORLD_SIZE/2, WORLD_SIZE/2); 00062 m_worldLV = new G4LogicalVolume (m_worldSol, air, "World", 0, 0, 0); 00063 m_worldPV = new G4PVPlacement (0, G4ThreeVector (), "World", m_worldLV, 0, false, 0); 00064 00065 // - Target 00066 G4ThreeVector posTarget (0, 0, -(TARGET_SIZE + TRACKER_SIZE)/2); 00067 m_targetSol = new G4Box ("target", TARGET_SIZE/2, TARGET_SIZE/2, TARGET_SIZE/2); 00068 m_targetLV = new G4LogicalVolume (m_targetSol, m_matPb, "Target",0,0,0); 00069 m_targetPV = new G4PVPlacement (0, posTarget, "Target", m_targetLV, m_worldPV, false, 0); 00070 00071 // - Tracker 00072 G4ThreeVector posTracker = G4ThreeVector (0,0,0); 00073 m_trackerSol = new G4Box ("tracker", TRACKER_SIZE/2, TRACKER_SIZE/2, TRACKER_SIZE/2); 00074 m_trackerLV = new G4LogicalVolume (m_trackerSol, air, "Tracker",0,0,0); 00075 m_trackerPV = new G4PVPlacement (0, posTracker, "Tracker", m_trackerLV, m_worldPV, false, 0); 00076 00077 // - Tracker segments: An example of Parameterised volumes; dummy 00078 // values for G4Box -- modified by parameterised volume 00079 m_chamberSol = new G4Box ("chamber", 100*cm, 100*cm, 10*cm); 00080 m_chamberLV = new G4LogicalVolume (m_chamberSol, m_matXenon, "Chamber",0,0,0); 00081 m_chamberPV = new G4PVParameterised ("Chamber", m_chamberLV, m_trackerPV, 00082 kZAxis, N_CHAMBERS, new VisG4ExampleParametrisation 00083 (N_CHAMBERS, -TRACKER_SIZE/2 + CHAMBER_SIZE/2, 00084 CHAMBER_SPACING, CHAMBER_SIZE, 00085 TRACKER_SIZE/10, TRACKER_SIZE)); 00086 00087 // - Sensitive detectors 00088 G4SDManager *sdm = G4SDManager::GetSDMpointer (); 00089 VisG4ExampleSD *sd = new VisG4ExampleSD ("TrackerChamberSD"); 00090 sdm->AddNewDetector (sd); 00091 m_chamberLV->SetSensitiveDetector (sd); 00092 00093 // - Visualization attributes 00094 G4VisAttributes *white= new G4VisAttributes (G4Colour (1.0, 1.0, 1.0)); 00095 G4VisAttributes *yellow = new G4VisAttributes (G4Colour (1.0, 1.0, 0.0)); 00096 00097 m_worldLV->SetVisAttributes (white); 00098 m_targetLV->SetVisAttributes (white); 00099 m_trackerLV->SetVisAttributes (white); 00100 m_chamberLV->SetVisAttributes (yellow); 00101 00102 return m_worldPV; 00103 }
double VisG4ExampleDetector::fullLength | ( | void | ) | const [inline] |
Definition at line 72 of file VisG4ExampleDetector.h.
References WORLD_SIZE.
Referenced by VisG4ExampleGenerator::GeneratePrimaries().
00073 { return WORLD_SIZE; }
const double VisG4ExampleDetector::CHAMBER_SIZE = 20 * cm [static, private] |
const double VisG4ExampleDetector::CHAMBER_SPACING = 80 * cm [static, private] |
G4LogicalVolume* VisG4ExampleDetector::m_chamberLV [private] |
G4VPhysicalVolume* VisG4ExampleDetector::m_chamberPV [private] |
G4Box* VisG4ExampleDetector::m_chamberSol [private] |
VisG4ExampleField* VisG4ExampleDetector::m_field [private] |
G4Material* VisG4ExampleDetector::m_matPb [private] |
G4Material* VisG4ExampleDetector::m_matXenon [private] |
G4LogicalVolume* VisG4ExampleDetector::m_targetLV [private] |
G4VPhysicalVolume* VisG4ExampleDetector::m_targetPV [private] |
G4Box* VisG4ExampleDetector::m_targetSol [private] |
G4LogicalVolume* VisG4ExampleDetector::m_trackerLV [private] |
G4VPhysicalVolume* VisG4ExampleDetector::m_trackerPV [private] |
G4Box* VisG4ExampleDetector::m_trackerSol [private] |
G4LogicalVolume* VisG4ExampleDetector::m_worldLV [private] |
G4VPhysicalVolume* VisG4ExampleDetector::m_worldPV [private] |
G4Box* VisG4ExampleDetector::m_worldSol [private] |
const int VisG4ExampleDetector::N_CHAMBERS = 5 [static, private] |
const double VisG4ExampleDetector::TARGET_SIZE = 5 * cm [static, private] |
const double VisG4ExampleDetector::TRACKER_SIZE = (N_CHAMBERS+1) * CHAMBER_SPACING [static, private] |
const double VisG4ExampleDetector::WORLD_SIZE = 1.2 * (TARGET_SIZE + TRACKER_SIZE) [static, private] |