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

#include <ExGflashDetectorConstruction.hh>

Inheritance diagram for ExGflashDetectorConstruction:
G4VUserDetectorConstruction

Public Member Functions

 ExGflashDetectorConstruction ()
 
 ~ExGflashDetectorConstruction ()
 
virtual G4VPhysicalVolumeConstruct ()
 
virtual void ConstructSDandField ()
 
const G4VPhysicalVolumeGetCristal (int num__crystal)
 
- 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

Definition at line 47 of file ExGflashDetectorConstruction.hh.

Constructor & Destructor Documentation

ExGflashDetectorConstruction::ExGflashDetectorConstruction ( )

Definition at line 70 of file ExGflashDetectorConstruction.cc.

References python.hepunit::cm, G4cout, and G4endl.

71  :fExperimentalHall_log(0),
72  fCalo_log(0),
73  fExperimentalHall_phys(0),
74  fCalo_phys(0),
75  fTheParameterisation(0),
76  fTheHMaker(0),
77  fTheParticleBounds(0),
78  fTheFastShowerModel(0)
79 {
80  G4cout<<"ExGflashDetectorConstruction::Detector constructor"<<G4endl;
81 
82  // Simplified `CMS-like` PbWO4 crystal calorimeter
83  fCalo_xside=31*cm;
84  fCalo_yside=31*cm;
85  fCalo_zside=24*cm;
86 
87  // GlashStuff
88  //Energy Cuts to kill particles
89 // fTheParticleBounds = new GFlashParticleBounds();
90  // Makes the EnergieSpots
91 // fTheHMaker = new GFlashHitMaker();
92 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
ExGflashDetectorConstruction::~ExGflashDetectorConstruction ( )

Definition at line 96 of file ExGflashDetectorConstruction.cc.

97 {
98  //@@@ ExGflashDetectorConstruction::Soll ich alles dlete
99 
100  // -- !! this is not properly deleted in MT where
101  // -- !! as there is one parameterisation/thread
102  // -- !! and only the last one is remembered.
103  if ( fTheParameterisation ) delete fTheParameterisation;
104  if ( fTheParticleBounds ) delete fTheParticleBounds;
105  if ( fTheHMaker ) delete fTheHMaker;
106  if ( fTheFastShowerModel ) delete fTheFastShowerModel;
107 }

Member Function Documentation

G4VPhysicalVolume * ExGflashDetectorConstruction::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Definition at line 111 of file ExGflashDetectorConstruction.cc.

References G4Region::AddRootLogicalVolume(), python.hepunit::cm, G4NistManager::FindOrBuildMaterial(), G4cout, G4endl, G4NistManager::Instance(), G4VisAttributes::Invisible, n, G4LogicalVolume::SetRegion(), and G4LogicalVolume::SetVisAttributes().

112 {
113  //--------- Definitions of Solids, Logical Volumes, Physical Volumes ---------
114  G4cout << "Defining the materials" << G4endl;
115  // Get nist material manager
116  G4NistManager* nistManager = G4NistManager::Instance();
117  // Build materials
118  G4Material* air = nistManager->FindOrBuildMaterial("G4_AIR");
119  G4Material* pbWO4 = nistManager->FindOrBuildMaterial("G4_PbWO4");
120 
121 
122  /*******************************
123  * The Experimental Hall *
124  *******************************/
125  fExperimentalHall_x=1000.*cm;
126  fExperimentalHall_y=1000.*cm;
127  fExperimentalHall_z=1000.*cm;
128 
129  fExperimentalHall_box = new G4Box("expHall_box", // World Volume
130  fExperimentalHall_x, // x size
131  fExperimentalHall_y, // y size
132  fExperimentalHall_z); // z size
133 
134  fExperimentalHall_log = new G4LogicalVolume(fExperimentalHall_box,
135  air,
136  "expHall_log",
137  0, //opt: fieldManager
138  0, //opt: SensitiveDetector
139  0); //opt: UserLimits
140  fExperimentalHall_phys = new G4PVPlacement(0,
141  G4ThreeVector(), //at (0,0,0)
142  "expHall",
143  fExperimentalHall_log,
144  0,
145  false,
146  0);
147 
148 
149  //------------------------------
150  // Calorimeter segments
151  //------------------------------
152  // Simplified `CMS-like` PbWO4 crystal calorimeter
153 
154  fNbOfCrystals = 10; // this are the crystals PER ROW in this example
155  // cube of 10 x 10 crystals
156  // don't change it @the moment, since
157  // the readout in event action assumes this
158  // dimensions and is not automatically adapted
159  // in this version of the example :-(
160  fCrystalWidth = 3*cm;
161  fCrystalLenght= 24*cm;
162  fCalo_xside=(fCrystalWidth*fNbOfCrystals)+1*cm;
163  fCalo_yside=(fCrystalWidth*fNbOfCrystals)+1*cm;
164  fCalo_zside=fCrystalLenght;
165 
166  G4Box *calo_box= new G4Box("CMS calorimeter", // its name
167  fCalo_xside/2., // size
168  fCalo_yside/2.,
169  fCalo_zside/2.);
170  fCalo_log = new G4LogicalVolume(calo_box, // its solid
171  air, // its material
172  "calo log", // its name
173  0, // opt: fieldManager
174  0, // opt: SensitiveDetector
175  0); // opt: UserLimit
176 
177  G4double Xpos = 0.0;
178  G4double Ypos = 0.0;
179  G4double Zpos = 100.0*cm;
180 
181  fCalo_phys = new G4PVPlacement(0,
182  G4ThreeVector(Xpos,Ypos,Zpos),
183  fCalo_log,
184  "calorimeter",
185  fExperimentalHall_log,
186  false,
187  1);
188  //Visibility
189  for (int i=0; i<fNbOfCrystals;i++)
190  {
191 
192  for (int j=0; j<fNbOfCrystals;j++)
193  {
194  int n = i*10+j;
195  fCrystal[n]= new G4Box("Crystal", // its name
196  fCrystalWidth/2,
197  fCrystalWidth/2,
198  fCrystalLenght/2); // size
199  fCrystal_log[n] = new G4LogicalVolume(fCrystal[n], // its solid
200  pbWO4, // its material
201  "Crystal_log"); // its name
202  G4ThreeVector crystalPos((i*fCrystalWidth)-135,
203  (j*fCrystalWidth)-135,0 );
204  fCrystal_phys[n] = new G4PVPlacement(0, // no rotation
205  crystalPos, // translation
206  fCrystal_log[n],
207  "crystal", // its name
208  fCalo_log,
209  false,
210  1);
211  }
212  }
213  G4cout << "There are " << fNbOfCrystals <<
214  " crystals per row in the calorimeter, so in total "<<
215  fNbOfCrystals*fNbOfCrystals << " crystals" << G4endl;
216  G4cout << "The have widthof " << fCrystalWidth /cm <<
217  " cm and a lenght of " << fCrystalLenght /cm
218  <<" cm. The Material is "<< pbWO4 << G4endl;
219 
220 
221  fExperimentalHall_log->SetVisAttributes(G4VisAttributes::Invisible);
222  G4VisAttributes* CaloVisAtt = new G4VisAttributes(G4Colour(1.0,1.0,1.0));
223  G4VisAttributes* CrystalVisAtt = new G4VisAttributes(G4Colour(1.0,1.0,0.0));
224  fCalo_log->SetVisAttributes(CaloVisAtt);
225  for (int i=0; i<100;i++)
226  {
227  fCrystal_log[i]->SetVisAttributes(CrystalVisAtt);
228  }
229  // define the parameterisation region
230  fRegion = new G4Region("crystals");
231  fCalo_log->SetRegion(fRegion);
232  fRegion->AddRootLogicalVolume(fCalo_log);
233 
234  return fExperimentalHall_phys;
235 }
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
CLHEP::Hep3Vector G4ThreeVector
void AddRootLogicalVolume(G4LogicalVolume *lv)
Definition: G4Region.cc:254
Definition: G4Box.hh:63
static G4NistManager * Instance()
void SetRegion(G4Region *reg)
G4GLOB_DLL std::ostream G4cout
const G4int n
static const G4VisAttributes Invisible
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetVisAttributes(const G4VisAttributes *pVA)
void ExGflashDetectorConstruction::ConstructSDandField ( )
virtual

Reimplemented from G4VUserDetectorConstruction.

Definition at line 239 of file ExGflashDetectorConstruction.cc.

References G4SDManager::AddNewDetector(), G4NistManager::FindOrBuildMaterial(), G4cout, G4endl, G4SDManager::GetSDMpointer(), G4NistManager::Instance(), GFlashShowerModel::SetHitMaker(), GFlashShowerModel::SetParameterisation(), GFlashShowerModel::SetParticleBounds(), and G4LogicalVolume::SetSensitiveDetector().

240 {
241  // -- sensitive detectors:
244  new ExGflashSensitiveDetector("Calorimeter",this);
245  SDman->AddNewDetector(CaloSD);
246 
247  for (int i=0; i<100;i++)
248  {
249  fCrystal_log[i]->SetSensitiveDetector(CaloSD);
250  }
251 
252  // Get nist material manager
253  G4NistManager* nistManager = G4NistManager::Instance();
254  G4Material* pbWO4 = nistManager->FindOrBuildMaterial("G4_PbWO4");
255  // -- fast simulation models:
256  // **********************************************
257  // * Initializing shower modell
258  // ***********************************************
259  G4cout << "Creating shower parameterization models" << G4endl;
260  GFlashShowerModel* lTheFastShowerModel =
261  new GFlashShowerModel("fastShowerModel",fRegion);
262  GFlashHomoShowerParameterisation* lTheParameterisation =
264  lTheFastShowerModel->SetParameterisation(*lTheParameterisation);
265  // Energy Cuts to kill particles:
266  GFlashParticleBounds* lTheParticleBounds = new GFlashParticleBounds();
267  lTheFastShowerModel->SetParticleBounds(*lTheParticleBounds);
268  // Makes the EnergieSpots
269  GFlashHitMaker* lTheHMaker = new GFlashHitMaker();
270  lTheFastShowerModel->SetHitMaker(*lTheHMaker);
271  G4cout<<"end shower parameterization."<<G4endl;
272  // **********************************************
273 }
void SetHitMaker(GFlashHitMaker &Maker)
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
void SetParticleBounds(GFlashParticleBounds &SpecificBound)
static G4NistManager * Instance()
G4GLOB_DLL std::ostream G4cout
void AddNewDetector(G4VSensitiveDetector *aSD)
Definition: G4SDManager.cc:67
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
#define G4endl
Definition: G4ios.hh:61
void SetParameterisation(GVFlashShowerParameterisation &DP)
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
const G4VPhysicalVolume* ExGflashDetectorConstruction::GetCristal ( int  num__crystal)
inline

Definition at line 57 of file ExGflashDetectorConstruction.hh.

Referenced by ExGflashSensitiveDetector::ProcessHits().

58  {return fCrystal_phys[num__crystal];};

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