FairRoot
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
FairModule.h
Go to the documentation of this file.
1 /********************************************************************************
2  * Copyright (C) 2014 GSI Helmholtzzentrum fuer Schwerionenforschung GmbH *
3  * *
4  * This software is distributed under the terms of the *
5  * GNU Lesser General Public Licence (LGPL) version 3, *
6  * copied verbatim in the file "LICENSE" *
7  ********************************************************************************/
8 #ifndef FAIRMODULE_H
9 #define FAIRMODULE_H
10 
11 #include "FairGeoInterface.h" // for FairGeoInterface
12 #include "FairGeoLoader.h" // for FairGeoLoader
13 #include "FairGeoNode.h" // for FairGeoNode
14 #include "FairGeoVolume.h" // for FairGeoVolume
15 #include "FairLogger.h"
16 #include "FairRun.h" // for FairRun
17 #include "FairRuntimeDb.h" // for FairRuntimeDb
18 
19 #include <Rtypes.h> // for Bool_t, Int_t, etc
20 #include <TList.h> // for TList (ptr only), TListIter
21 #include <TNamed.h> // for TNamed
22 #include <TObjArray.h> // for TObjArray
23 #include <TString.h> // for TString, operator!=
24 #include <string> // for string
25 
26 class FairVolumeList;
27 class FairVolume;
28 class TArrayI;
29 class TGeoMatrix;
30 class TGeoNode;
31 class TGeoVolume;
32 class TRefArray;
33 class TVirtualMC;
34 
46 class FairModule : public TNamed
47 {
48  public:
50  FairModule();
52  FairModule(const char* Name, const char* title, Bool_t Active = kFALSE);
54  virtual ~FairModule();
56  virtual void Print(Option_t*) const { ; }
58  virtual void SetGeometryFileName(TString fname, TString geoVer = "0");
60  virtual TString GetGeometryFileName() { return fgeoName; }
62  virtual TString GetGeometryFileVer() { return fgeoVer; }
64  virtual void ConstructGeometry();
66  virtual void ConstructOpGeometry();
68  virtual void ConstructRootGeometry(TGeoMatrix* shiftM = nullptr);
70  virtual void ConstructASCIIGeometry();
72  virtual void ModifyGeometry()
73  __attribute__((deprecated("Use FairAlignmentHandler instead, see Tutorial4 for examples")))
74  {
75  LOG(warn) << "This function is deprecated. Use FairAlignmentHandler instead, see Tutorial4 for examples.";
76  }
77  virtual void RegisterAlignmentMatrices() { ; }
78 
80  virtual void ConstructGDMLGeometry(__attribute__((unused)) TGeoMatrix* posrot);
83  virtual void SetSpecialPhysicsCuts() { ; }
85  virtual FairModule* CloneModule() const;
87  virtual void BeginWorkerRun() const { ; }
89  virtual void FinishWorkerRun() const { ; }
90 
92  template<class T, class U>
93  void ConstructASCIIGeometry(T* dataType1, TString containerName = "", U* datatype2 = nullptr);
94 
97  virtual Bool_t IsSensitive(const std::string& name);
99  virtual Bool_t CheckIfSensitive(__attribute__((unused)) std::string name) __attribute__((
100  deprecated("The method CheckIfSensitive is deprecated. Implement IsSensitive in the detector classes.")))
101  {
102  return kFALSE;
103  }
105  virtual void ExpandNode(TGeoNode* Node);
107  virtual void ExpandNodeForGDML(__attribute__((unused)) TGeoNode* curNode);
109  virtual Int_t getVolId(const TString&) const { return 0; }
111  Int_t GetModId() { return fModId; }
113  void SetVerboseLevel(Int_t level) { fVerboseLevel = level; }
115  Bool_t IsActive() { return fActive; }
117  void SetModId(Int_t id) { fModId = id; }
121  void SetMotherVolume(TString volName) { fMotherVolumeName = volName; }
123  void ProcessNodes(TList* aList);
125  virtual void SetParContainers() { ; }
128  virtual void InitParContainers() { ; }
130  TList* GetListOfGeoPar() { return flGeoPar; }
131 
133  static thread_local FairVolumeList* vList;
134 
135  static thread_local Int_t fNbOfVolumes;
136 
137  static thread_local TRefArray* svList;
138 
139  static thread_local TArrayI* volNumber;
142  void AddSensitiveVolume(TGeoVolume* v);
143 
144  private:
146  void SetDefaultMatrixName(TGeoMatrix* matrix);
147  void AssignMediumAtImport(TGeoVolume* v); // O.Merle, 29.02.2012 - see impl.
148 
150  void ReAssignMediaId();
151  void swap(FairModule& other) throw();
152 
153  protected:
154  FairModule(const FairModule&);
156  TString fgeoVer;
157  TString fgeoName;
158  Int_t fModId;
159  Bool_t fActive;
162  TList* flGeoPar;
163  Bool_t fGeoSaved;
164  TVirtualMC* fMC;
165 
166  ClassDef(FairModule, 4);
167 };
168 
169 template<class T, class U>
170 void FairModule::ConstructASCIIGeometry(T* dataType1, TString containerName, U*)
171 {
173  FairGeoInterface* GeoInterface = loader->getGeoInterface();
174  T* MGeo = new T();
175  MGeo->print();
176  MGeo->setGeomFile(GetGeometryFileName());
177  GeoInterface->addGeoModule(MGeo);
178  Bool_t rc = GeoInterface->readSet(MGeo);
179  if (rc) {
180  MGeo->create(loader->getGeoBuilder());
181  }
182 
183  TList* volList = MGeo->getListOfVolumes();
184  // store geo parameter
185  FairRun* fRun = FairRun::Instance();
187 
188  dataType1 = MGeo;
189 
190  if ("" != containerName) {
191  LOG(info) << "Add GeoNodes for " << MGeo->getDescription() << " to container " << containerName;
192 
193  // U par=(U)(rtdb->getContainer(containerName));
194  U* par = static_cast<U*>(rtdb->getContainer(containerName));
195  TObjArray* fSensNodes = par->GetGeoSensitiveNodes();
196  TObjArray* fPassNodes = par->GetGeoPassiveNodes();
197 
198  TListIter iter(volList);
199  FairGeoNode* node = nullptr;
200  FairGeoVolume* aVol = nullptr;
201 
202  while ((node = static_cast<FairGeoNode*>(iter.Next()))) {
203  aVol = dynamic_cast<FairGeoVolume*>(node);
204  if (node->isSensitive()) {
205  fSensNodes->AddLast(aVol);
206  } else {
207  fPassNodes->AddLast(aVol);
208  }
209  }
210  ProcessNodes(volList);
211  par->setChanged();
212  par->setInputVersion(fRun->GetRunId(), 1);
213  }
214 }
215 
216 #endif // FAIRMODULE_H
virtual void SetGeometryFileName(TString fname, TString geoVer="0")
Definition: FairModule.cxx:199
void AddSensitiveVolume(TGeoVolume *v)
Definition: FairModule.cxx:294
virtual void see Tutorial4 for examples LOG(warn)<< "This function is deprecated. Use FairAlignmentHandler instead
virtual void ConstructGeometry()
Definition: FairModule.cxx:60
list of container factories
Definition: FairRuntimeDb.h:24
for(Int_t i=0;i< 100;i++)
Bool_t readSet(FairGeoSet *)
static thread_local TRefArray * svList
Definition: FairModule.h:137
virtual void ConstructOpGeometry()
Definition: FairModule.cxx:66
static thread_local Int_t fNbOfVolumes
Definition: FairModule.h:135
void SetModId(Int_t id)
Definition: FairModule.h:117
virtual void RegisterAlignmentMatrices()
Definition: FairModule.h:77
virtual Int_t getVolId(const TString &) const
Definition: FairModule.h:109
Bool_t isSensitive()
Definition: FairGeoNode.h:146
virtual void SetParContainers()
Definition: FairModule.h:125
virtual void SetSpecialPhysicsCuts()
Definition: FairModule.h:83
void ProcessNodes(TList *aList)
Definition: FairModule.cxx:244
FairGeoInterface * getGeoInterface()
Definition: FairGeoLoader.h:34
static FairGeoLoader * Instance()
static FairRun * Instance()
Definition: FairRun.cxx:31
Int_t fVerboseLevel
Definition: FairModule.h:161
void addGeoModule(FairGeoSet *)
virtual ~FairModule()
Definition: FairModule.cxx:72
Bool_t fGeoSaved
list of Detector Geometry parameters
Definition: FairModule.h:163
virtual void FinishWorkerRun() const
Definition: FairModule.h:89
Bool_t fActive
Definition: FairModule.h:159
Int_t fModId
Definition: FairModule.h:158
virtual void Print(Option_t *) const
Definition: FairModule.h:56
virtual void ConstructGDMLGeometry(__attribute__((unused)) TGeoMatrix *posrot)
Definition: FairModule.cxx:508
static thread_local FairVolumeList * vList
Definition: FairModule.h:133
virtual void ConstructASCIIGeometry()
Definition: FairModule.cxx:552
TVirtualMC * fMC
flag for initialisation
Definition: FairModule.h:164
virtual void ExpandNode(TGeoNode *Node)
Definition: FairModule.cxx:569
static thread_local TArrayI * volNumber
Definition: FairModule.h:139
FairParSet * getContainer(const Text_t *)
virtual void ExpandNodeForGDML(__attribute__((unused)) TGeoNode *curNode)
Definition: FairModule.cxx:516
FairModule & operator=(const FairModule &)
Definition: FairModule.cxx:149
TList * GetListOfGeoPar()
Definition: FairModule.h:130
ClassDef(FairModule, 4)
cahed pointer to MC (available only after initialization)
FairMQExParamsParOne * par
Int_t fNbOfSensitiveVol
Definition: FairModule.h:160
TList * flGeoPar
Definition: FairModule.h:162
virtual FairModule * CloneModule() const
Definition: FairModule.cxx:711
virtual TString GetGeometryFileVer()
Definition: FairModule.h:62
void SetMotherVolume(TString volName)
Definition: FairModule.h:121
FairGeoBuilder * getGeoBuilder()
Definition: FairGeoLoader.h:35
virtual void ConstructRootGeometry(TGeoMatrix *shiftM=nullptr)
Definition: FairModule.cxx:327
virtual void InitParContainers()
Definition: FairModule.h:128
FairRuntimeDb * GetRuntimeDb(void)
Definition: FairRun.h:80
virtual TString GetGeometryFileName()
Definition: FairModule.h:60
TString fgeoName
Definition: FairModule.h:157
virtual void see Tutorial4 for examples see Tutorial4 for examples
Definition: FairModule.h:75
virtual Bool_t CheckIfSensitive(__attribute__((unused)) std::string name) __attribute__((deprecated("The method CheckIfSensitive is deprecated. Implement IsSensitive in the detector classes.")))
Definition: FairModule.h:99
void SetVerboseLevel(Int_t level)
Definition: FairModule.h:113
virtual void ModifyGeometry() __attribute__((deprecated("Use FairAlignmentHandler instead
virtual void BeginWorkerRun() const
Definition: FairModule.h:87
virtual Bool_t IsSensitive(const std::string &name)
Definition: FairModule.cxx:559
TString fMotherVolumeName
Definition: FairModule.h:140
FairVolume * getFairVolume(FairGeoNode *fNode)
Definition: FairModule.cxx:312
Bool_t IsActive()
Definition: FairModule.h:115
Int_t GetRunId()
Definition: FairRun.h:97
Int_t GetModId()
Definition: FairModule.h:111
TString fgeoVer
Definition: FairModule.h:156