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