]> git.uio.no Git - u/mrichter/AliRoot.git/blob - AliGeant4/AliSingleModuleConstruction.cxx
added THijing, CONTAINERS categories
[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   fVersion = right.fVersion;
49   fType = right.fType;
50   fProcessConfig = right.fProcessConfig;
51   fAllLVSensitive = right.fAllLVSensitive;
52   fSDManager = right.fSDManager;
53
54   G4String moduleName = right.fModuleName;
55   moduleName.toLower();
56   fMessenger = new AliSingleModuleConstructionMessenger(this, moduleName);
57 }
58
59 AliSingleModuleConstruction::AliSingleModuleConstruction() {
60 //
61 }
62
63 AliSingleModuleConstruction::~AliSingleModuleConstruction() {
64 //
65   delete fMessenger;
66 }
67
68 // operators
69
70 AliSingleModuleConstruction& 
71 AliSingleModuleConstruction::operator=(const AliSingleModuleConstruction& right)
72 {    
73   // check assignement to self
74   if (this == &right) return *this;
75   
76   // base class assignement
77   AliModuleConstruction::operator=(right);
78   
79   fVersion = right.fVersion;
80   fType = right.fType;
81   fProcessConfig = right.fProcessConfig;
82   fAllLVSensitive = right.fAllLVSensitive;
83   fSDManager = right.fSDManager;
84   
85   return *this;
86 }
87
88 // private methods
89
90 void AliSingleModuleConstruction::CreateSensitiveDetectors()
91 {
92 // Creates sensitive detectors.
93 // ---
94
95   if (fAllLVSensitive)
96     CreateSensitiveDetectors1();
97   else
98     CreateSensitiveDetectors2();
99
100   // set static number of logical volumes already processed
101   G4LogicalVolumeStore* pLVStore = G4LogicalVolumeStore::GetInstance();
102   fSDManager->SetNofLVWithSD(pLVStore->entries());  
103 }   
104
105 void AliSingleModuleConstruction::CreateSensitiveDetectors1()
106
107 // Creates sensitive detectors.
108 // Sensitive detectors are set to all logical volumes
109 // ---
110
111   G4LogicalVolumeStore* pLVStore = G4LogicalVolumeStore::GetInstance();
112   G4int nofLV = pLVStore->entries();
113   
114   G4int nofLVWithSD = fSDManager->GetNofLVWithSD();
115   
116   for (G4int i=nofLVWithSD; i<nofLV; i++) {
117     G4LogicalVolume* lv = (*pLVStore)[i];
118     fSDManager->CreateSD(lv, fAliModule);
119   }
120 }
121
122 void AliSingleModuleConstruction::CreateSensitiveDetectors2()
123
124 // Creates sensitive detectors.
125 // Sensitive detectors are set only to logical volumes
126 // in G3SensVolVector.
127 // ---
128
129   TG4GeometryManager* pGeometryManager = TG4GeometryManager::Instance();
130
131   G3SensVolVector pSVVector
132     = pGeometryManager->GetG3SensVolVector();
133
134   G4int nofSV = pSVVector.entries();
135   if (nofSV>0)
136     for (G4int isv=0; isv<nofSV; isv++) {
137       G4LogicalVolume* lv = pSVVector[isv];
138       fSDManager->CreateSD(lv, fAliModule);
139     } 
140 }
141
142 // public methods 
143
144 void AliSingleModuleConstruction::Configure()
145
146 // Executes the detector setup Root macro
147 // (extracted from AliRoot Config.C) and
148 // G4 macro.
149 // ---
150           
151   // filepaths and macro names 
152   G4String rootFilePath;
153   G4String g4FilePath;
154
155   if (fType == kDetector) {
156     // macro file path 
157     //$AG4_INSTALL/macro/XXX/Config.C
158     //$AG4_INSTALL/macro/XXX/Config.in
159     rootFilePath = AliFiles::DetConfig1();
160     rootFilePath.append(fModuleName);
161     rootFilePath.append(AliFiles::DetConfig2());
162     g4FilePath = rootFilePath;
163
164     rootFilePath.append(AliFiles::DetConfig3());
165     g4FilePath.append(AliFiles::DetConfig4());
166
167     // path to g3calls.dat file
168     //$AG4_INSTALL/macro/XXX/g3calls_vN.dat
169     fDataFilePath = AliFiles::DetData1();
170     fDataFilePath.append(fModuleName);
171     fDataFilePath.append(AliFiles::DetData2());
172     G4String version("v");
173     AliGlobals::AppendNumberToString(version,fVersion);
174     fDataFilePath.append(version);
175     fDataFilePath.append(AliFiles::DetData3());           
176   }  
177   else if (fType == kStructure) {    
178     // macro file path is set to 
179     //$AG4_INSTALL/macro/STRUCT/XXXConfig.C
180     //$AG4_INSTALL/macro/STRUCT/XXXConfig.in
181     rootFilePath = AliFiles::DetConfig1();
182     rootFilePath.append(AliFiles::STRUCT());
183     rootFilePath.append(AliFiles::DetConfig2());
184     rootFilePath.append(fModuleName);
185     g4FilePath = rootFilePath;
186
187     rootFilePath.append(AliFiles::DetConfig3());
188     g4FilePath.append(AliFiles::DetConfig4());
189
190     // path to g3calls.dat file
191     //$AG4_INSTALL/macro/STRUCT/g3calls_XXXvN.dat
192     fDataFilePath = AliFiles::DetData1();
193     fDataFilePath.append(AliFiles::STRUCT());
194     fDataFilePath.append(AliFiles::DetData2());
195     fDataFilePath.append(fModuleName);
196     G4String version("v");
197     AliGlobals::AppendNumberToString(version,fVersion);
198     fDataFilePath.append(version);
199     fDataFilePath.append(AliFiles::DetData3());
200   }  
201   
202   if (fProcessConfig) {
203     // load and execute aliroot macro
204     gROOT->LoadMacro(rootFilePath);
205     G4String macroName = AliFiles::DetConfigName1();
206     AliGlobals::AppendNumberToString(macroName, fVersion);
207     macroName.append(AliFiles::DetConfigName2());
208     gInterpreter->ProcessLine(macroName);
209   } 
210   
211   // add process g4 config macro
212   // execute Geant4 macro if file is specified as an argument 
213   G4String command = "/control/execute ";
214   G4UImanager* pUI = G4UImanager::GetUIpointer();  
215   pUI->ApplyCommand(command + g4FilePath);
216   
217   // get AliModule created in Config.C macro
218   fAliModule = gAlice->GetModule(fModuleName);
219   if (!fAliModule) {
220     G4String text = "AliSingleModuleConstruction::Configure:\n";
221     text = text + "    AliModule " + fModuleName;
222     text = text + " has not been found in gAlice.";
223     AliGlobals::Exception(text);
224   }  
225 }
226
227 void AliSingleModuleConstruction::Construct()
228
229 // Constructs geometry.
230 // ---
231
232   // print default element table
233   // const G4ElementTable* table = G4Element::GetElementTable();
234   // G4cout << "Default elemnt table: " << G4endl;
235   // for (G4int i=0; i<table->entries(); i++) {
236   //   G4cout << *(*table)[i] << G4endl;
237   // }  
238
239   // Configure();
240
241   // get geometry manager
242   TG4GeometryManager* pGeometryManager = TG4GeometryManager::Instance();
243
244   // register module name in the name map
245   pGeometryManager->SetMapSecond(fAliModule->GetName());        
246
247   if (fReadGeometry) {
248     // create G3 geometry from g3calls.dat
249     pGeometryManager->SetWriteGeometry(false);
250     pGeometryManager->ReadG3Geometry(fDataFilePath);
251   }
252   else {
253     // set geometry output stream for this module
254     pGeometryManager->SetWriteGeometry(fWriteGeometry);
255     if (fWriteGeometry) 
256       pGeometryManager->OpenOutFile(fDataFilePath);
257
258     // create geometry from AliRoot
259
260     // construct materials
261     fAliModule->CreateMaterials();
262
263     // construct G3 geometry
264     fAliModule->CreateGeometry();
265         
266     if (fWriteGeometry) 
267       pGeometryManager->CloseOutFile();
268   }  
269   
270   // construct G4 geometry
271   G4VPhysicalVolume* world = pGeometryManager->CreateG4Geometry();
272   if (!fgWorld) fgWorld = world; 
273   
274   // set the detector frame (envelope)
275   // (without warning output if enevelope is not defined)
276   SetDetFrame(false);
277
278   // create sensitive detectors
279   CreateSensitiveDetectors();
280
281   // build sensitive detectors table
282   fAliModule->Init();
283
284   // construct geometry for display
285   fAliModule->BuildGeometry();
286
287   // reset TG4GeometryManager 
288   pGeometryManager->ClearG3Tables();
289   
290   // print current total number of logical volumes
291   G4cout << "Current total number of sensitive volumes: "
292          << pGeometryManager->NofVolumes() << G4endl;
293
294 #ifdef ALICE_VISUALIZE
295   if (GetDetFrame()) {
296     // set visualization attributes
297     // if detector envelope is defined
298     SetDetVisibility(true);
299     SetDetColour("Yellow");
300   }  
301 #endif
302 }