KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVGeoImport.cpp
Go to the documentation of this file.
1 //Created by KVClassFactory on Tue Apr 23 13:04:21 2013
2 //Author: John Frankland,,,
3 
4 #include "KVGeoImport.h"
5 #include <KVMultiDetArray.h>
7 #include <TGeoBBox.h>
8 #include <TPluginManager.h>
9 
10 #include <KVEvent.h>
11 #include <KVGroup.h>
12 #include <TGeoPhysicalNode.h>
14 #include <KVNamedParameter.h>
15 #include <KVDataAnalyser.h>
16 
18 
19 
20 
21 
29 KVGeoImport::KVGeoImport(TGeoManager* g, KVIonRangeTable* r, KVMultiDetArray* m, Bool_t create) : KVGeoNavigator(g), fCreateArray(create), fOrigin(nullptr)
30 {
31  // Import geometry described by TGeoManager into KVMultiDetArray object
32  //
33  // if create=kFALSE, we do not create any detectors etc., but just set up
34  // the required links between the geometry and the existing array detectors
35  //
36  // We change the colours of volumes depending on their material
37 
38  fArray = m;
39  fRangeTable = r;
40  fCurrentTrajectory.SetAddToNodes(kFALSE);
41  fCurrentTrajectory.SetPathInTitle(kFALSE);
42  fCheckDetVolNames = kFALSE;
43 
44  std::map<std::string, EColor> colors;
45  colors["Si"] = EColor(kGray + 1);
46  colors["CsI"] = EColor(kBlue - 8);
47  colors["Al"] = EColor(kGray);
48  colors["Cu"] = EColor(kOrange + 2);
49  TGeoVolume* vol;
50  TIter next(g->GetListOfVolumes());
51  while ((vol = (TGeoVolume*)next())) {
52  TGeoMedium* med = vol->GetMedium();
53  if (!med) continue;
54  TGeoMaterial* mat = med->GetMaterial();
55  if (colors.find(mat->GetName()) != colors.end()) vol->SetLineColor(colors[mat->GetName()]);
56  if (mat->GetDensity() < 0.1) vol->SetTransparency(60);
57  }
58 }
59 
60 
61 
64 
66 {
67  // Destructor
69 }
70 
71 
72 
75 
77 {
78  // All detectors crossed by the particle's trajectory are added to the multidetector
79 
80  KVDetector* detector = GetCurrentDetector();
81  if (!detector) return;
82  detector->GetNode()->SetName(detector->GetName());
83  if (fLastDetector && detector != fLastDetector) {
84  fLastDetector->GetNode()->AddBehind(detector);
85  detector->GetNode()->AddInFront(fLastDetector);
87  }
88  fLastDetector = detector;
89 }
90 
91 
92 
98 
100  Double_t ThetaMin, Double_t PhiMin, Double_t ThetaMax, Double_t PhiMax)
101 {
102  // Scan the geometry in order to find all detectors and detector alignments.
103  // This is done by sending out "particles" from (0,0,0) or (x,y,z) (if SetOrigin(x,y,z) was called)
104  // in all directions between (ThetaMin,ThetaMax) - with respect to Z-axis - and (PhiMin,PhiMax) - cylindrical
105  // angle in the (X,Y)-plane, over a grid of step dTheta in Theta and dPhi in Phi.
106 
107  Int_t ndets0 = fArray->GetDetectors()->GetEntries();
108 
109  unique_ptr<KVEvent> evt(new KVEvent());
110  KVNucleus* nuc = evt->AddParticle();
111  nuc->SetZAandE(1, 1, 1);
112 
113  Info("ImportGeometry",
114  "Importing geometry in angular ranges : Theta=[%f,%f:%f] Phi=[%f,%f:%f]", ThetaMin, ThetaMax, dTheta, PhiMin, PhiMax, dPhi);
115  if (fOrigin)
116  Info("ImportGeometry",
117  "Origin for geometry = (%f,%f,%f)", fOrigin->x(), fOrigin->y(), fOrigin->z());
118  Int_t count = 0;
120  std::cout << "\xd" << "Info in <KVGeoImport::ImportGeometry>: tested " << count << " directions" << std::flush;
121  for (Double_t theta = ThetaMin; theta <= ThetaMax; theta += dTheta) {
122  for (Double_t phi = PhiMin; phi <= PhiMax; phi += dPhi) {
123  nuc->SetTheta(theta);
124  nuc->SetPhi(phi);
125  fLastDetector = nullptr;
126  PropagateEvent(evt.get(), fOrigin);
127  count++;
129  std::cout << "\xd" << "Info in <KVGeoImport::ImportGeometry>: tested " << count << " directions" << std::flush;
130  }
131  }
133  std::cout << "Info in <KVGeoImport::ImportGeometry>: tested " << count << " directions" << std::endl;
134  else
135  std::cout << std::endl;
136 
137  Info("ImportGeometry",
138  "Imported %d detectors into array", fArray->GetDetectors()->GetEntries() - ndets0);
139 
140  // make sure detector nodes are correct
141  TIter next(fArray->GetDetectors());
142  KVDetector* d;
143  while ((d = (KVDetector*)next())) {
144  d->GetNode()->RehashLists();
145  d->SetNameOfArray(fArray->GetName());
146  }
147  // set up all detector node trajectories
149  // create all groups
151 
152  if (fCreateArray) {
154  // Set up internal navigator of array with all informations on detector/
155  // structure name formatting, correspondance lists, etc.
158  for (int i = 0; i < fStrucNameFmt.GetEntries(); i++) {
160  nav->SetStructureNameFormat(fmt->GetName(), fmt->GetString());
161  }
163  nav->AbsorbDetectorPaths(this);
167  }
168 }
169 
170 
171 
173 
175 {
176  fLastDetector = d;
177 }
178 
179 
180 
184 
186 {
187  // Override KVGeoNavigator method
188  // We build the list of all trajectories through the array
189 
191 
192  KVGeoNavigator::PropagateParticle(nuc, TheOrigin);
193 
194  if (fLastDetector && fLastDetector->GetNode()) {
195  if (fCurrentTrajectory.GetN()) {
198  }
200  }
201  else {
203  }
206  fArray->AddTrajectory(tr);
207  }
208  }
209 }
210 
211 
212 
217 
219 {
220  // Add to list of acceptable detector names when scanning geometry
221  // Each detector volume name will be tested; if it doesn't contain
222  // any of the (partial) detector names in the list, it will be ignored
223 
226 }
227 
228 
229 
236 
238 {
239  // Returns pointer to KVDetector corresponding to current location
240  // in geometry. Detector is created and added to array if needed.
241  // We also set up any geometry structure elements (from nodes beginning with "STRUCT_")
242  // Any detector volume with a name not recognised by comparison with the list
243  // fAcceptedDetectorNames (if defined) will be ignored.
244 
245  KVString detector_name;
246  Bool_t multilay;
247  TGeoVolume* detector_volume = GetCurrentDetectorNameAndVolume(detector_name, multilay);
248  // failed to identify current volume as part of a detector
249 
250  if (!detector_volume) return 0;
251  // check volume belongs to us
252  if (fCheckDetVolNames) {
253  int nn = fAcceptedDetectorNames.GetNpar();
254  bool reject_volume = true;
255  for (int i = 0; i < nn; ++i) {
256  if (detector_name.Contains(fAcceptedDetectorNames.GetNameAt(i))) {
257  reject_volume = false;
258  break;
259  }
260  }
261  if (reject_volume) return 0;
262  }
263 
264  // has detector already been built ? if not, do it now
265  KVDetector* det = fArray->GetDetector(detector_name);
266  if (!fCreateArray) {
267  if (det) {
268  // set matrix & shape for entrance window if not done yet
269  if (!det->GetEntranceWindow().ROOTGeo()) {
271  det->SetEntranceWindowShape((TGeoBBox*)GetCurrentVolume()->GetShape());
272  }
273  TString vol_name(GetCurrentVolume()->GetName());
274  if (!multilay || vol_name.BeginsWith("ACTIVE_")) {
275  // set matrix & shape for active layer
277  det->SetActiveLayerShape((TGeoBBox*)GetCurrentVolume()->GetShape());
278  // set full path to physical node as title of detector's node (KVGeoDetectorNode)
279  det->GetNode()->SetTitle(GetCurrentPath());
280  }
281  // add entry to correspondance list between physical nodes and detectors (all layers)
283  }
284  }
285  else {
286  if (!det) {
287  det = BuildDetector(detector_name, detector_volume);
288  if (det) {
289  // Setting the entrance window shape and matrix
290  // ============================================
291  // for consistency, the matrix and shape MUST correspond
292  // i.e. we cannot have the matrix corresponding to the entrance window
293  // of a multilayer detector and the shape corresponding to the
294  // whole detector (all layers) - otherwise, calculation of points
295  // on detector entrance window will be false!
296 // Info("GetCurrentDetector","Setting EW matrix to current matrix:");
297 // GetCurrentMatrix()->Print();
299  det->SetEntranceWindowShape((TGeoBBox*)GetCurrentVolume()->GetShape());
300  TString vol_name(GetCurrentVolume()->GetName());
301  if (!multilay || vol_name.BeginsWith("ACTIVE_")) {
302  // first layer of detector (or only layer) is also active layer
303 // Info("GetCurrentDetector","and also setting active layer matrix to current matrix:");
304 // GetCurrentMatrix()->Print();
306  det->SetActiveLayerShape((TGeoBBox*)GetCurrentVolume()->GetShape());
307  // set full path to physical node as title of detector's node (KVGeoDetectorNode)
308  det->GetNode()->SetTitle(GetCurrentPath());
309  }
310  // add entry to correspondance list between physical nodes and detectors (all layers)
312  fArray->Add(det);
313  Int_t nstruc = CurrentStructures().GetEntries();
314  if (nstruc) {
315  // Build and add geometry structure elements
316  KVGeoStrucElement* ELEM = fArray;
317  for (int i = 0; i < nstruc; i++) {
319  KVGeoStrucElement* nextELEM = ELEM->GetStructure(elem->GetName());
320  if (!nextELEM) {
321  // make new structure
322  nextELEM = new KVGeoStrucElement(elem->GetName(), elem->GetType());
323  nextELEM->SetNumber(elem->GetNumber());
324  ELEM->Add(nextELEM);
325  }
326  ELEM = nextELEM;
327  }
328  // add detector to last structure
329  ELEM->Add(det);
330  }
331  }
332  }
333  else {
334  // Detector already built, are we now in its active layer ?
335  TString vol_name(GetCurrentVolume()->GetName());
336  if (!multilay || vol_name.BeginsWith("ACTIVE_")) {
337 // Info("GetCurrentDetector","Setting active layer matrix to current matrix:");
338 // GetCurrentMatrix()->Print();
340  det->SetActiveLayerShape((TGeoBBox*)GetCurrentVolume()->GetShape());
341  // set full path to physical node as title of detector's node (KVGeoDetectorNode)
342  det->GetNode()->SetTitle(GetCurrentPath());
343  }
344  // add entry to correspondance list between physical nodes and detectors (all layers)
346  }
347  }
348  return det;
349 }
350 
351 
352 
379 
381 {
382  // Create a KVDetector with given name for the given volume
383  //
384  // #### Detector definition in geometry ####
385  // 1. All detector volumes & nodes must have names which begin with **DET_**
386  //
387  // 2. They must be made of materials which are known by the range table fRangeTable.
388  //
389  // 3. For multi-layer detectors, the "active" layer volume/node must have a name beginning with **ACTIVE_**
390  //
391  // 4. The "thickness" of the detector or any layer inside a multilayer detector will be taken as the size of the volume's
392  // shape along its Z-axis (so make sure that you define your detector volumes in this way).
393  //
394  // 5. It is assumed that the natural length units of the geometry are centimetres.
395  //
396  // 6. The name of the KVDetector object created and added to the array will be taken
397  // from the unique full path of the node corresponding to the geometrical positioning
398  // of the detector, see KVGeoNavigator::ExtractDetectorNameFromPath()
399  //
400  // 7. The 'type' of the detector will be set to the name of the material
401  // in the detector's active layer i.e. if active layer material name is "Si",
402  // detector type will be 'Si'
403  //
404  // 8. Default class for all detectors is KVDetector. If you want to use another class
405  // you need to define it using SetDetectorPlugin() method and put the
406  // associated line in your .kvrootrc configuration file.
407 
408 
409  KVDetector* d = 0;
410  TPluginHandler* ph = NULL;
411  if (fDetectorPlugin == "" || !(ph = LoadPlugin("KVDetector", fDetectorPlugin))) {
412  d = new KVDetector;
413  }
414  else {
415  d = (KVDetector*)ph->ExecPlugin(0);
416  }
417 
418  d->SetName(det_name);
419 
420  Int_t nlayer = det_vol->GetNdaughters();
421  if (nlayer) {
422  for (int i = 0; i < nlayer; i++) {
423  AddLayer(d, det_vol->GetNode(i)->GetVolume());
424  }
425  }
426  else
427  AddLayer(d, det_vol);
428  TString type = d->GetActiveLayer()->GetName();
429  //type.ToUpper();
430  d->SetType(type);
431  return d;
432 }
433 
434 
435 
441 
443 {
444  // Add an absorber layer to the detector.
445  //
446  // Volumes representing 'active' layers in detectors must have names
447  // which begin with **ACTIVE_**.
448 
449  TString vnom = vol->GetName();
450  // exclude dead zone layers
451  if (vnom.BeginsWith("DEADZONE")) return;
452  TGeoMaterial* material = vol->GetMaterial();
453  KVIonRangeTableMaterial* irmat = fRangeTable->GetMaterial(material);
454  if (!irmat) {
455  Warning("AddLayer", "Unknown material %s/%s used in layer %s of detector %s",
456  material->GetName(), material->GetTitle(), vol->GetName(), det->GetName());
457  return;
458  }
459  TGeoBBox* sh = dynamic_cast<TGeoBBox*>(vol->GetShape());
460  if (!sh) {
461  Warning("AddLayer", "Unknown shape class %s used in layer %s of detector %s",
462  vol->GetShape()->ClassName(), vol->GetName(), det->GetName());
463  return; // just in case - for now, all shapes derive from TGeoBBox...
464  }
465  Double_t width = 2.*sh->GetDZ(); // thickness in centimetres
466  KVMaterial* absorber;
467  if (irmat->IsGas()) {
468  Double_t p = material->GetPressure();
469  Double_t T = material->GetTemperature();
470  absorber = new KVMaterial(irmat->GetType(), width, p, T);
471  }
472  else
473  absorber = new KVMaterial(irmat->GetType(), width);
474  det->AddAbsorber(absorber);
475  if (vnom.BeginsWith("ACTIVE_")) det->SetActiveLayer(absorber);
476 }
477 
478 
479 
int Int_t
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
ROOT::R::TRInterface & r
#define SafeDelete(p)
#define d(i)
const Bool_t kFALSE
bool Bool_t
double Double_t
const Bool_t kTRUE
EColor
kGray
kOrange
kBlue
include TDocParser_001 C image html pict1_TDocParser_001 png width
int type
Color * colors
virtual void SetNumber(UInt_t num)
Definition: KVBase.h:209
const Char_t * GetType() const
Definition: KVBase.h:170
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition: KVBase.cpp:756
UInt_t GetNumber() const
Definition: KVBase.h:213
static Bool_t IsRunningBatchAnalysis()
Base class for detector geometry description, interface to energy-loss calculations.
Definition: KVDetector.h:121
void AddAbsorber(KVMaterial *)
Definition: KVDetector.cpp:753
const KVPosition & GetEntranceWindow() const
Definition: KVDetector.h:647
void SetActiveLayer(KVMaterial *actif)
Definition: KVDetector.h:244
void SetActiveLayerMatrix(const TGeoHMatrix *)
KVGeoDetectorNode * GetNode()
Definition: KVDetector.h:285
void SetEntranceWindowMatrix(const TGeoHMatrix *)
Set ROOT geometry global matrix transformation to coordinate frame of entrance window.
void SetEntranceWindowShape(TGeoBBox *)
Set ROOT geometry shape of entrance window.
void SetActiveLayerShape(TGeoBBox *)
Set ROOT geometry shape of active layer volume.
Base class container for multi-particle events.
Definition: KVEvent.h:176
Path taken by particles through multidetector geometry.
Bool_t Contains(const Char_t *name) const
void ReverseOrder()
Reverse the order of the nodes in the trajectory.
void AddLast(KVGeoDetectorNode *n)
void Clear(Option_t *="")
Clear list of nodes in trajectory.
void AddInFront(KVDetector *)
void AddBehind(KVDetector *)
Import detector array described by ROOT geometry and set up corresponding KVMultiDetArray object.
Definition: KVGeoImport.h:67
KVNameValueList fAcceptedDetectorNames
Definition: KVGeoImport.h:79
KVDetector * GetCurrentDetector()
virtual void ParticleEntersNewVolume(KVNucleus *)
All detectors crossed by the particle's trajectory are added to the multidetector.
Definition: KVGeoImport.cpp:76
KVDetector * BuildDetector(TString det_name, TGeoVolume *det_vol)
KVGeoDNTrajectory fCurrentTrajectory
Definition: KVGeoImport.h:73
Bool_t fCheckDetVolNames
Definition: KVGeoImport.h:80
void AddAcceptedDetectorName(const char *name)
TString fDetectorPlugin
Definition: KVGeoImport.h:72
KVDetector * fLastDetector
Definition: KVGeoImport.h:70
void PropagateParticle(KVNucleus *, TVector3 *TheOrigin=nullptr)
virtual ~KVGeoImport()
Destructor.
Definition: KVGeoImport.cpp:65
KVIonRangeTable * fRangeTable
Definition: KVGeoImport.h:69
void AddLayer(KVDetector *, TGeoVolume *)
KVMultiDetArray * fArray
Definition: KVGeoImport.h:68
TVector3 * fOrigin
Definition: KVGeoImport.h:74
void SetLastDetector(KVDetector *)
Bool_t fCreateArray
Definition: KVGeoImport.h:71
void ImportGeometry(Double_t dTheta=0.1, Double_t dPhi=1.0, Double_t ThetaMin=0.0, Double_t PhiMin=0.0, Double_t ThetaMax=180.0, Double_t PhiMax=360.0)
Definition: KVGeoImport.cpp:99
Link physical path to node in geometry with detector.
Base class for propagation of particles through array geometry.
void SetDetectorNameFormat(const Char_t *fmt)
void PropagateEvent(KVEvent *, TVector3 *TheOrigin=0)
const TGeoHMatrix * GetCurrentMatrix() const
Returns pointer to internal copy of current global transformation matrix.
void AbsorbDetectorPaths(KVGeoNavigator *GN)
KVUniqueNameList fDetectorPaths
correspondance between physical node and detector objects
KVNameValueList fStrucNameFmt
list of user-defined formats for structure names
const TClonesArray & CurrentStructures() const
void SetNameCorrespondanceList(const Char_t *)
TEnv * fDetStrucNameCorrespList
list(s) of correspondance for renaming structures/detectors
virtual void PropagateParticle(KVNucleus *, TVector3 *TheOrigin=0)
TString GetCurrentPath() const
TGeoVolume * GetCurrentVolume() const
KVString fDetNameFmt
user-defined format for detector names
void SetStructureNameFormat(const Char_t *type, const Char_t *fmt)
TGeoVolume * GetCurrentDetectorNameAndVolume(KVString &, Bool_t &)
TGeoManager * GetGeometry() const
Base class describing elements of array geometry.
KVDetector * GetDetector(const Char_t *name) const
Return detector in this structure with given name.
const KVSeqCollection * GetDetectors() const
virtual void Add(KVBase *)
KVGeoStrucElement * GetStructure(const Char_t *name) const
Material for use in energy loss & range calculations.
Abstract base class for calculation of range & energy loss of charged particles in matter.
KVIonRangeTableMaterial * GetMaterial(const Char_t *material) const
Returns pointer to material of given name or type.
Description of physical materials used to construct detectors; interface to range tables.
Definition: KVMaterial.h:41
Base class for describing the geometry of a detector array.
void DeduceIdentificationTelescopesFromGeometry()
void CalculateDetectorSegmentationIndex()
void AddTrajectory(KVGeoDNTrajectory *d)
KVGeoNavigator * GetNavigator() const
void SetGeometry(TGeoManager *)
const TSeqCollection * GetTrajectories() const
void CalculateReconstructionTrajectories()
void AssociateTrajectoriesAndNodes()
void DeduceGroupsFromTrajectories()
KVNamedParameter * GetParameter(Int_t idx) const
return the parameter object with index idx
void SetValue(const Char_t *name, value_type value)
const Char_t * GetNameAt(Int_t idx) const
Int_t GetNpar() const
return the number of stored parameters
Int_t GetEntries() const
A generic named parameter storing values of different types.
const Char_t * GetString() const
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
void SetZAandE(Int_t z, Int_t a, Double_t ekin)
Set atomic number, mass number, and kinetic energy in MeV.
Definition: KVNucleus.cpp:733
void SetTheta(Double_t theta)
Definition: KVParticle.h:654
void SetPhi(Double_t phi)
Definition: KVParticle.h:658
Bool_t ROOTGeo() const
Returns kTRUE if ROOT geometry is used, kFALSE if not.
Definition: KVPosition.cpp:598
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
virtual void Add(TObject *obj)
virtual TObject * FindObject(const char *name) const
virtual Int_t GetEntries() const
virtual Double_t GetDZ() const
Double_t GetTemperature() const
Double_t GetPressure() const
virtual Double_t GetDensity() const
TGeoMaterial * GetMaterial() const
TGeoVolume * GetVolume() const
TGeoShape * GetShape() const
TGeoMedium * GetMedium() const
Int_t GetNdaughters() const
void SetTransparency(Char_t transparency=0)
TGeoNode * GetNode(const char *name) const
virtual void SetLineColor(Color_t lcolor)
TGeoMaterial * GetMaterial() const
virtual const char * GetName() const
virtual void SetTitle(const char *title="")
virtual const char * GetTitle() const
virtual void SetName(const char *name)
Int_t GetEntries() const
virtual const char * ClassName() const
virtual void Warning(const char *method, const char *msgfmt,...) const
virtual void Info(const char *method, const char *msgfmt,...) const
Longptr_t ExecPlugin(int nargs, const T &... params)
Bool_t BeginsWith(const char *s, ECaseCompare cmp=kExact) const
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Double_t z() const
Double_t x() const
Double_t y() const
const long double m
Definition: KVUnits.h:70
const long double g
masses
Definition: KVUnits.h:72
double T(double x)