]> git.uio.no Git - u/mrichter/AliRoot.git/blob - AliGeant4/AliSingleModuleConstruction.cxx
private method GetSDConstruction() added; AliSDManager usage replaced with AliSDConst...
[u/mrichter/AliRoot.git] / AliGeant4 / AliSingleModuleConstruction.cxx
1 // $Id$
2 // Category: geometry
3 //
4 // See the class description in the header file.
5
6 #include "AliSingleModuleConstruction.h"
7 #include "AliSingleModuleConstructionMessenger.h"
8 #include "AliSDManager.h"
9 #include "AliGlobals.h"
10 #include "AliFiles.h"
11 #include "AliRun.h"
12 #include "AliModule.h"
13
14 #include "TG4GeometryManager.h"
15
16 #include <G3SensVolVector.hh>
17 #include <G4UImanager.hh>
18 //#include <G4Element.hh>
19 #include <G4LogicalVolume.hh>
20 #include <G4LogicalVolumeStore.hh>
21
22 #include <TROOT.h> 
23 #include <TCint.h> 
24
25 G4VPhysicalVolume* AliSingleModuleConstruction::fgWorld = 0;
26
27 AliSingleModuleConstruction::AliSingleModuleConstruction(
28                                 G4String moduleName, G4int version,
29                                 AliModuleType moduleType)
30   : AliModuleConstruction(moduleName),
31     fVersion(version),
32     fType(moduleType),
33     fProcessConfig(true),
34     fAllLVSensitive(false)
35 {
36 //
37   fSDManager = AliSDManager::Instance();
38
39   moduleName.toLower();
40   fMessenger = new AliSingleModuleConstructionMessenger(this, moduleName);
41 }
42
43 AliSingleModuleConstruction::AliSingleModuleConstruction(
44                                 const AliSingleModuleConstruction& right)
45   : AliModuleConstruction(right)                                
46 {
47 //
48   // allocation in assignement operator
49   fMessenger = 0;
50   
51   // copy stuff
52   *this = right;
53 }
54
55 AliSingleModuleConstruction::AliSingleModuleConstruction() {
56 //
57 }
58
59 AliSingleModuleConstruction::~AliSingleModuleConstruction() {
60 //
61   delete fMessenger;
62 }
63
64 // operators
65
66 AliSingleModuleConstruction& 
67 AliSingleModuleConstruction::operator=(const AliSingleModuleConstruction& right)
68 {    
69   // check assignement to self
70   if (this == &right) return *this;
71   
72   // base class assignement
73   AliModuleConstruction::operator=(right);
74   
75   fVersion = right.fVersion;
76   fType = right.fType;
77   fProcessConfig = right.fProcessConfig;
78   fAllLVSensitive = right.fAllLVSensitive;
79   fSDManager = right.fSDManager;
80  
81   // messenger is protected from copying  
82   if (fMessenger) delete fMessenger;
83   G4String moduleName = right.fModuleName;
84   moduleName.toLower();
85   fMessenger = new AliSingleModuleConstructionMessenger(this, moduleName);
86  
87   return *this;
88 }
89
90 // private methods
91
92 void AliSingleModuleConstruction::CreateSensitiveDetectors()
93 {
94 // Creates sensitive detectors.
95 // ---
96
97   if (fAllLVSensitive)
98     CreateSensitiveDetectors1();
99   else
100     CreateSensitiveDetectors2();
101
102   // set static number of logical volumes already processed
103   G4LogicalVolumeStore* pLVStore = G4LogicalVolumeStore::GetInstance();
104   fSDManager->SetNofLVWithSD(pLVStore->size());  
105 }   
106
107 void AliSingleModuleConstruction::CreateSensitiveDetectors1()
108
109 // Creates sensitive detectors.
110 // Sensitive detectors are set to all logical volumes
111 // ---
112
113   G4LogicalVolumeStore* pLVStore = G4LogicalVolumeStore::GetInstance();
114   G4int nofLV = pLVStore->size();
115   
116   G4int nofLVWithSD = fSDManager->GetNofLVWithSD();
117   
118   for (G4int i=nofLVWithSD; i<nofLV; i++) {
119     G4LogicalVolume* lv = (*pLVStore)[i];
120     fSDManager->CreateSD(lv, fAliModule);
121   }
122 }
123
124 void AliSingleModuleConstruction::CreateSensitiveDetectors2()
125
126 // Creates sensitive detectors.
127 // Sensitive detectors are set only to logical volumes
128 // in G3SensVolVector.
129 // ---
130
131   TG4GeometryManager* pGeometryManager = TG4GeometryManager::Instance();
132
133   G3SensVolVector pSVVector
134     = pGeometryManager->GetG3SensVolVector();
135
136   G4int nofSV = pSVVector.entries();
137   if (nofSV>0)
138     for (G4int isv=0; isv<nofSV; isv++) {
139       G4LogicalVolume* lv = pSVVector[isv];
140       fSDManager->CreateSD(lv, fAliModule);
141     } 
142 }
143
144 // public methods 
145
146 void AliSingleModuleConstruction::Configure(const AliFiles& files)
147
148 // Executes the detector setup Root macro
149 // (extracted from AliRoot Config.C) and
150 // G4 macro.
151 // ---
152           
153   // filepaths and macro names 
154   G4bool isStructure = (fType == kStructure);
155   G4String rootFilePath 
156     = files.GetRootMacroPath(fModuleName, isStructure);
157   G4String g4FilePath
158     = files.GetG4MacroPath(fModuleName, isStructure);
159   fDataFilePath 
160     = files.GetG3CallsDatPath(fModuleName, fVersion, isStructure); 
161   
162   // load and execute aliroot config macro
163   if (fProcessConfig) {
164     gROOT->LoadMacro(rootFilePath);
165     G4String macroName = files.GetDefaultMacroName();
166     //macroName = macroName + "_" + fModuleName;
167     macroName = macroName + "(";
168     AliGlobals::AppendNumberToString(macroName, fVersion);
169     macroName = macroName + ")";
170     gInterpreter->ProcessLine(macroName);
171   } 
172   
173   // process g4 config macro
174   G4String command = "/control/execute ";
175   G4UImanager* pUI = G4UImanager::GetUIpointer();  
176   pUI->ApplyCommand(command + g4FilePath);
177   
178   // get AliModule created in Config.C macro
179   fAliModule = gAlice->GetModule(fModuleName);
180   if (!fAliModule) {
181     G4String text = "AliSingleModuleConstruction::Configure:\n";
182     text = text + "    AliModule " + fModuleName;
183     text = text + " has not been found in gAlice.";
184     AliGlobals::Exception(text);
185   }  
186 }
187
188 void AliSingleModuleConstruction::Construct()
189
190 // Constructs geometry.
191 // ---
192
193   // print default element table
194   // const G4ElementTable* table = G4Element::GetElementTable();
195   // G4cout << "Default elemnt table: " << G4endl;
196   // for (G4int i=0; i<table->entries(); i++) {
197   //   G4cout << *(*table)[i] << G4endl;
198   // }  
199
200   // Configure();
201
202   // get geometry manager
203   TG4GeometryManager* pGeometryManager = TG4GeometryManager::Instance();
204
205   // register module name in the name map
206   pGeometryManager->SetMapSecond(fAliModule->GetName());        
207
208   if (fReadGeometry) {
209     // create G3 geometry from g3calls.dat
210     pGeometryManager->SetWriteGeometry(false);
211     pGeometryManager->ReadG3Geometry(fDataFilePath);
212   }
213   else {
214     // set geometry output stream for this module
215     pGeometryManager->SetWriteGeometry(fWriteGeometry);
216     if (fWriteGeometry) 
217       pGeometryManager->OpenOutFile(fDataFilePath);
218
219     // create geometry from AliRoot
220
221     // construct materials
222     fAliModule->CreateMaterials();
223
224     // construct G3 geometry
225     fAliModule->CreateGeometry();
226         
227     if (fWriteGeometry) 
228       pGeometryManager->CloseOutFile();
229   }  
230   
231   // construct G4 geometry
232   G4VPhysicalVolume* world = pGeometryManager->CreateG4Geometry();
233   if (!fgWorld) fgWorld = world; 
234   
235   // set the detector frame (envelope)
236   // (without warning output if enevelope is not defined)
237   SetDetFrame(false);
238
239   // create sensitive detectors
240   CreateSensitiveDetectors();
241
242   // build sensitive detectors table
243   fAliModule->Init();
244
245   // construct geometry for display
246   fAliModule->BuildGeometry();
247
248   // reset TG4GeometryManager 
249   pGeometryManager->ClearG3Tables();
250   
251   // print current total number of logical volumes
252   G4cout << "Current total number of sensitive volumes: "
253          << pGeometryManager->NofVolumes() << G4endl;
254
255 #ifdef ALICE_VISUALIZE
256   if (GetDetFrame()) {
257     // set visualization attributes
258     // if detector envelope is defined
259     SetDetVisibility(true);
260     SetDetColour("Yellow");
261   }  
262 #endif
263 }