]> git.uio.no Git - u/mrichter/AliRoot.git/blob - AliGeant4/AliModulesComposition.cxx
updated commands description
[u/mrichter/AliRoot.git] / AliGeant4 / AliModulesComposition.cxx
1 // $Id$
2 // Category: geometry
3 //
4 // Author: I. Hrivnacova
5 //
6 // Class AliModulesComposition
7 // ---------------------------
8 // See the class description in the header file.
9
10 #include "AliModulesComposition.h"
11 #include "AliSingleModuleConstruction.h"
12 #include "AliMoreModulesConstruction.h"
13 #include "AliDetSwitch.h"
14 #include "AliMagneticField.h"
15 #include "AliGlobals.h"
16 #include "AliFiles.h"
17
18 #include "TG4XMLGeometryGenerator.h"
19 #include "TG4GeometryManager.h"
20
21 #include <G4Material.hh>
22 #include <G4VPhysicalVolume.hh>
23
24 //_____________________________________________________________________________
25 AliModulesComposition::AliModulesComposition()
26   : fReadGeometry(false),
27     fWriteGeometry(false),
28     fMagneticField(0),
29     fMessenger(this)
30 {
31 //
32   fMoreModulesConstruction = new AliMoreModulesConstruction();
33 }
34
35 //_____________________________________________________________________________
36 AliModulesComposition::AliModulesComposition(const AliModulesComposition& right)
37   : fMessenger(this)
38 {
39 //
40   AliGlobals::Exception("AliModulesComposition is protected from copying.");  
41 }
42
43 //_____________________________________________________________________________
44 AliModulesComposition::~AliModulesComposition() {
45 //   
46   delete fMoreModulesConstruction;
47   delete fMagneticField;
48   
49   // destroy det switch vector
50   DetSwitchIterator it;
51   for (it = fDetSwitchVector.begin(); it != fDetSwitchVector.end(); it++)
52     delete *it; 
53   
54   // destroy det construction vector
55   SingleModuleIterator itm;
56   for (itm = fModuleConstructionVector.begin(); 
57        itm != fModuleConstructionVector.end(); it++)
58     delete *itm;
59 }
60
61 // operators
62
63 //_____________________________________________________________________________
64 AliModulesComposition& 
65 AliModulesComposition::operator=(const AliModulesComposition& right)
66 {
67   // check assignement to self
68   if (this == &right) return *this;
69   
70   AliGlobals::Exception("AliModulesComposition is protected from assigning.");  
71
72   return *this;  
73 }    
74           
75 // protected methods
76
77 //_____________________________________________________________________________
78 void AliModulesComposition::AddDetSwitch(AliDetSwitch* detSwitch)
79 {
80 // Adds detSwitch to the detSwitch vector.
81 // ---
82
83   fDetSwitchVector.push_back(detSwitch);
84   fMessenger.SetCandidates();
85 }  
86   
87 //_____________________________________________________________________________
88 void AliModulesComposition::AddSingleModuleConstruction(
89                                const G4String& name, G4int version, 
90                                AliModuleType moduleType)
91 {
92 // Adds SingleModuleConstruction.
93 // ---
94
95   fModuleConstructionVector
96     .push_back(new AliSingleModuleConstruction(name, version, moduleType));
97 }  
98                   
99 //_____________________________________________________________________________
100 void AliModulesComposition::AddMoreModuleConstruction(
101                                const G4String& name, G4int version, 
102                                AliModuleType moduleType)
103 {
104 // Adds module to MoreModulesConstruction (construction of dependent
105 // modules.)
106 // ---
107
108   fMoreModulesConstruction->AddModule(name, version, moduleType);
109 }  
110                                   
111 //_____________________________________________________________________________
112 void AliModulesComposition::ConstructModules()
113 {
114 // Construct geometry of all modules (both standalone and dependent.)
115 // ---
116
117   // set common options
118   SetReadGeometryToModules(fReadGeometry);
119   SetWriteGeometryToModules(fWriteGeometry);
120   
121   // configure single modules
122   SingleModuleIterator it;
123   for (it  = fModuleConstructionVector.begin(); 
124        it != fModuleConstructionVector.end(); it++) {
125        
126     (*it)->Configure(*AliFiles::Instance());
127     cout << "Module " << (*it)->GetDetName() << " configured." << endl;
128   }  
129      
130   // configure dependent modules
131   if (fMoreModulesConstruction->GetNofModules() > 0)
132     fMoreModulesConstruction->Configure(*AliFiles::Instance());
133
134   // construct single modules
135   for (it  = fModuleConstructionVector.begin(); 
136        it != fModuleConstructionVector.end(); it++) {
137
138     G4cout << "Module " << (*it)->GetDetName()
139            << " will be constructed now." << G4endl;
140     (*it)->Construct();
141   }  
142     
143   // construct dependent modules
144   if (fMoreModulesConstruction->GetNofModules() > 0) {
145     G4cout << "Dependent modules will be constructed now." << G4endl;
146     fMoreModulesConstruction->Construct();
147   }  
148 }  
149
150 //_____________________________________________________________________________
151 AliDetSwitch* AliModulesComposition::GetDetSwitch(
152                                         const G4String& moduleName) const
153 {
154 // Returns the detector switch with given detector name.
155 // ---
156
157   DetSwitchConstIterator it;
158   for (it = fDetSwitchVector.begin(); it != fDetSwitchVector.end(); it++)
159     if ((*it)->GetDetName() == moduleName) return *it; 
160
161   G4String text = "AliModulesComposition::GetDetSwitch:\n";
162   text = text + "Wrong detector name for " + moduleName;   
163   AliGlobals::Exception(text);
164   return 0;  
165
166
167 //_____________________________________________________________________________
168 void AliModulesComposition::SetReadGeometryToModules(G4bool readGeometry)
169 {
170 // Sets readGeometry control to all modules.
171 // ---
172
173   // single module constructions
174   SingleModuleIterator it;
175   for (it  = fModuleConstructionVector.begin(); 
176        it != fModuleConstructionVector.end(); it++)        
177     (*it)->SetReadGeometry(readGeometry);
178
179   // more modules construction
180   for (G4int i=0; i<fMoreModulesConstruction->GetNofModules(); i++) 
181     fMoreModulesConstruction
182       ->GetModuleConstruction(i)->SetReadGeometry(readGeometry);
183 }    
184   
185 //_____________________________________________________________________________
186 void AliModulesComposition::SetWriteGeometryToModules(G4bool writeGeometry)
187 {
188 // Sets writeGeometry control to all modules.
189 // ---
190
191   // single module constructions
192   SingleModuleIterator it;
193   for (it  = fModuleConstructionVector.begin(); 
194        it != fModuleConstructionVector.end(); it++)        
195     (*it)->SetWriteGeometry(writeGeometry);
196
197   // more modules construction
198   for (G4int i=0; i<fMoreModulesConstruction->GetNofModules(); i++) 
199     fMoreModulesConstruction
200       ->GetModuleConstruction(i)->SetWriteGeometry(writeGeometry);
201 }    
202
203 //_____________________________________________________________________________
204 void AliModulesComposition::SetProcessConfigToModules(G4bool processConfig)
205 {
206 // Sets processConfig control to all modules.
207 // ---
208
209   // single module constructions
210   SingleModuleIterator it;
211   for (it  = fModuleConstructionVector.begin(); 
212        it != fModuleConstructionVector.end(); it++)        
213     (*it)->SetProcessConfig(processConfig);
214   
215   // more modules construction
216   for (G4int i=0; i<fMoreModulesConstruction->GetNofModules(); i++) 
217     fMoreModulesConstruction
218       ->GetModuleConstruction(i)->SetProcessConfig(processConfig);
219 }    
220
221 // public methods
222
223 //_____________________________________________________________________________
224 void AliModulesComposition::SwitchDetOn(const G4String& moduleNameVer)
225
226 // Switchs on module specified by name and version.
227 // ---
228
229   DetSwitchIterator it;
230
231   if (moduleNameVer == "ALL") {
232     for (it = fDetSwitchVector.begin(); it != fDetSwitchVector.end(); it++)
233       (*it)->SwitchOnDefault(); 
234   }
235   else if (moduleNameVer == "NONE") {
236     for (it = fDetSwitchVector.begin(); it != fDetSwitchVector.end(); it++)
237       (*it)->SwitchOff(); 
238   }
239   else {
240     // get version number
241     G4int len = moduleNameVer.length();
242     G4String moduleName = moduleNameVer.substr(0, len-1);
243     G4String version = moduleNameVer.substr(len-1, 1);
244     G4int iVersion = AliGlobals::StringToInt(version);
245
246     if (iVersion < 0) {
247       // in case the version number is not provided
248       // the default one is set
249       SwitchDetOnDefault(moduleNameVer);
250     }  
251     else 
252       SwitchDetOn(moduleName, iVersion);
253   }
254 }
255
256 //_____________________________________________________________________________
257 void AliModulesComposition::SwitchDetOn(const G4String& moduleName, 
258                                         G4int version)
259
260 // Switchs on module specified by name and version.
261 // ---
262
263   GetDetSwitch(moduleName)->SwitchOn(version);
264 }
265
266 //_____________________________________________________________________________
267 void AliModulesComposition::SwitchDetOnDefault(const G4String& moduleName)
268
269 // Switchs on module specified by name with default version.
270 // ---
271
272   GetDetSwitch(moduleName)->SwitchOnDefault();
273 }
274
275 //_____________________________________________________________________________
276 void AliModulesComposition::SwitchDetOff(const G4String& moduleName)
277
278 // Switchs off module specified by name.
279 // ---
280
281   if (moduleName == "ALL") {
282     DetSwitchIterator it;
283     for (it = fDetSwitchVector.begin(); it != fDetSwitchVector.end(); it++)
284       (*it)->SwitchOff(); 
285   }
286   else 
287     GetDetSwitch(moduleName)->SwitchOff();
288 }
289
290 //_____________________________________________________________________________
291 void AliModulesComposition::PrintSwitchedDets() const
292
293 // Lists switched detectors.
294 // ---
295
296   G4String svList = GetSwitchedDetsList();
297     
298   G4cout << "Switched Alice detectors: " << G4endl;
299   G4cout << "--------------------------" << G4endl;
300   G4cout << svList << G4endl;
301 }
302
303 //_____________________________________________________________________________
304 void AliModulesComposition::PrintAvailableDets() const
305
306 // Lists available detectors.
307 // ---
308
309   G4String avList = GetAvailableDetsList();
310     
311   G4cout << "Available Alice detectors: " << G4endl;
312   G4cout << "---------------------------" << G4endl;
313   G4cout << avList << G4endl;
314 }
315
316 //_____________________________________________________________________________
317 void AliModulesComposition::PrintMaterials() const
318 {
319 // Prints all materials.
320 // ---
321
322   const G4MaterialTable* matTable = G4Material::GetMaterialTable();
323   G4cout << *matTable;
324 }
325
326 //_____________________________________________________________________________
327 void AliModulesComposition::GenerateXMLGeometry() const 
328 {
329 // Generates XML geometry file from the top volume.
330 // The file name is set according the last switched detector
331 // registered in the det switch vector.
332 // ---
333
334   G4VPhysicalVolume* world = AliSingleModuleConstruction::GetWorld();
335
336   // XML filename
337   // according to last switched detector
338   G4String detName;
339   G4String detVersion = "";
340   G4int version = -1;
341   for (G4int i=fDetSwitchVector.size()-1; i>=0; i--) {
342     version = fDetSwitchVector[i]->GetSwitchedVersion();
343     if (version > -1) {
344       detName = fDetSwitchVector[i]->GetDetName();
345       AliGlobals::AppendNumberToString(detVersion,version); 
346       break;
347     }  
348   }  
349   G4String filePath 
350     = AliFiles::Instance()->GetXMLFilePath(detName, version);
351   
352   // set top volume name
353   G4String topName = world->GetName() + "_comp";
354   
355   // generate XML
356   
357   TG4XMLGeometryGenerator xml;
358   xml.OpenFile(filePath);
359
360   // generate materials 
361   // not yet implemented
362   // xml.GenerateMaterials(version, "today", "Generated from G4",
363   //                     "v4", world->GetLogicalVolume());
364
365   // generate volumes tree
366   xml.GenerateSection(detName, detVersion, "today", "Generated from Geant4",
367                       topName, world->GetLogicalVolume());
368   xml.CloseFile();
369   
370   // set verbose
371   G4cout << "File " << detName << "v" << version << ".xml has been generated." 
372          << G4endl;
373 }  
374
375 //_____________________________________________________________________________
376 G4String AliModulesComposition::GetSwitchedDetsList() const
377
378 // Returns list of switched detectors.
379 // ---
380
381   G4String svList = "";  
382   G4int nofSwitchedDets = 0;
383   DetSwitchConstIterator it;
384   
385   for (it = fDetSwitchVector.begin(); it != fDetSwitchVector.end(); it++) {
386     G4int iVersion = (*it)->GetSwitchedVersion();
387     if (iVersion > -1) {
388       nofSwitchedDets++;
389       G4String moduleNameVer = (*it)->GetDetName();
390       AliGlobals::AppendNumberToString(moduleNameVer, iVersion);
391       svList += moduleNameVer;
392       svList += " "; 
393     }
394   }
395
396   if (nofSwitchedDets == fDetSwitchVector.size()) svList = "ALL: " + svList;
397   if (nofSwitchedDets == 0) svList = "NONE";   
398
399   return svList;
400 }
401
402 //_____________________________________________________________________________
403 G4String AliModulesComposition::GetAvailableDetsList() const
404
405 // Returns list of available detectors.
406 // ---
407
408   G4String svList = "";
409   DetSwitchConstIterator it;
410   
411   for (it = fDetSwitchVector.begin(); it != fDetSwitchVector.end(); it++)
412     for (G4int iv=0; iv<(*it)->GetNofVersions(); iv++) {
413       G4String moduleNameVer = (*it)->GetDetName();
414       AliGlobals::AppendNumberToString(moduleNameVer, iv);
415       svList += moduleNameVer;
416       svList += " ";
417     } 
418
419   return svList;
420 }
421
422 //_____________________________________________________________________________
423 G4String AliModulesComposition::GetAvailableDetsListWithCommas() const
424
425 // Returns list of available detectors with commas.
426 // ---
427
428   G4String svList = "";
429   G4int id =0;
430   DetSwitchConstIterator it;
431
432   for (it = fDetSwitchVector.begin(); it != fDetSwitchVector.end(); it++)
433     for (G4int iv=0; iv<(*it)->GetNofVersions(); iv++) {
434       G4String moduleNameVer = (*it)->GetDetName();
435       AliGlobals::AppendNumberToString(moduleNameVer, iv);
436       svList += moduleNameVer;
437       if (iv < (*it)->GetNofVersions()-1)        svList += "/";
438       else if (id++ < fDetSwitchVector.size()-1) svList += ", ";
439     }
440
441   return svList;
442 }
443
444 //_____________________________________________________________________________
445 G4String AliModulesComposition::GetDetNamesList() const
446
447 // Returns list of detector names.
448 // ---
449
450   G4String svList = "";
451   DetSwitchConstIterator it;
452   
453   for (it = fDetSwitchVector.begin(); it != fDetSwitchVector.end(); it++) {
454     svList += (*it)->GetDetName();
455     svList += " ";
456   }
457
458   return svList;
459 }
460
461 //_____________________________________________________________________________
462 G4String AliModulesComposition::GetDetNamesListWithCommas() const
463
464 // Returns list of detector names with commas.
465 // ---
466
467   G4String svList = "";
468   G4int id =0;
469   DetSwitchConstIterator it;
470
471   for (it = fDetSwitchVector.begin(); it != fDetSwitchVector.end(); it++) {
472     svList += (*it)->GetDetName();
473     if (id++ < fDetSwitchVector.size()-1) svList += ", ";
474   }
475
476   return svList;
477 }
478
479 //_____________________________________________________________________________
480 void AliModulesComposition::SetMagField(G4double fieldValue)
481 {
482 // Sets uniform magnetic field to specified value.
483 // ---
484
485   // create fields if it does not exist
486   if (!fMagneticField) fMagneticField = new AliMagneticField();
487   
488   // set value
489   fMagneticField->SetFieldValue(fieldValue);
490 }
491