]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVE/EveDet/AliEveFMDLoader.cxx
Merge of EVE-dev branch.
[u/mrichter/AliRoot.git] / EVE / EveDet / AliEveFMDLoader.cxx
CommitLineData
274799f3 1/**************************************************************************
2 * Copyright(c) 2008, Christian Holm Christensen *
3 * *
4 * Author: Christian Holm Christensen. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation for any purposes is hereby granted without fee, *
9 * provided that the above copyright notice appears in all copies and *
10 * that both the copyright notice and this permission notice remains *
11 * intact. The authors make no claims about the suitability of this *
12 * software for any purpose. It is provided "as is" without express or *
13 * implied warranty. *
14 **************************************************************************/
15/* $Id$ */
16/** @file AliEveFMDLoader.cxx
17 @author Christian Holm Christensen <cholm@nbi.dk>
18 @date Sun Mar 26 17:59:18 2006
19 @brief Implementation of AliEveFMDLoader singleton class
20*/
21//____________________________________________________________________
22//
23// Forward Multiplicity Detector based on Silicon wafers. This class
24// is the loader for the event display.
25//
26// This class is a singleton, meaning that there's only one instance
27// of this. This is done to speed up the processing by putting all
28// things that are needed every time into the constructor.
29//
30#include "AliRunLoader.h"
31#include "EveBase/AliEveEventManager.h"
32#include "AliEveFMDLoader.h"
33#include "../FMD/AliFMDUShortMap.h"
34#include "../FMD/AliFMDBoolMap.h"
35#include "../FMD/AliFMDGeometry.h"
36#include "../FMD/AliFMDParameters.h"
37#include "../FMD/AliFMDDetector.h"
38#include "../FMD/AliFMDRing.h"
39#include "../FMD/AliFMDBaseDigit.h"
40#include "../FMD/AliFMDDigit.h"
41#include "../FMD/AliFMDRawReader.h"
42#include "../FMD/AliFMDHit.h"
43#include "AliESDEvent.h"
44#include "AliESDFMD.h"
45#include "AliLog.h"
46#include "AliRawReader.h"
47#include <TClonesArray.h>
48#include <TTree.h>
49#include <TGeoShape.h>
50#include <TGeoManager.h>
51#include <TEveGeoNode.h>
52#include <TEveBoxSet.h>
53#include <TEveQuadSet.h>
54#include <TEveManager.h>
55#include <TEveUtil.h>
56#include <TStyle.h>
57#include <TMath.h>
58#include <iostream>
59
60// Some very private variables.
61namespace
62{
63 const Char_t* kDetector = "FMD%d";
64 const Char_t* kRing = "FMD%d%c";
65 const Char_t* kModule = "FMD%d%c[%02d-%02d]";
66 const Char_t* kSector = "FMD%d%c[%02d] %s";
67
68 const Char_t* kHits = "Hits";
69 const Char_t* kDigits = "Digits";
70 const Char_t* kRaw = "Raw";
71 const Char_t* kESD = "ESD";
72
73}
74
75//____________________________________________________________________
76AliEveFMDLoader* AliEveFMDLoader::fgInstance = 0;
77
78//____________________________________________________________________
79AliEveFMDLoader* AliEveFMDLoader::Instance()
80{
81 // Get the singleton instance. If the instance has not been
82 // instantised yet, it will be after this call.
83 if (!fgInstance)
84 fgInstance = new AliEveFMDLoader();
85 return fgInstance;
86}
87
88//____________________________________________________________________
89AliEveFMDLoader::AliEveFMDLoader(const char* name, Bool_t useBoxes,
90 Bool_t old)
91 : TEveElementList(name, 0),
92 fHitPalette(0, 1000),
93 fDigitPalette(0, 1023),
94 fMultPalette(0, 20),
95 fUseBoxDigits(useBoxes),
96 fHitCache("AliFMDHit",0),
97 fDigitCache("AliFMDDigit", 0),
98 fRawCache("AliFMDDigit", 0)
99{
100 // Constructor
101 // Parameters:
102 // @param name Name of the folder.
103 // @param useBoxes Whether to use boxes or Quads for the signals
104 // @param old Whether to enable reading old RCU data format
105
106 // increase reference count
107 IncDenyDestroy();
108
109 // Increase reference counts on palettes
110 fHitPalette.IncRefCount();
111 fDigitPalette.IncRefCount();
112 fMultPalette.IncRefCount();
113
114 // Initialize the FMD geometry manager
72c0b0d5 115 AliEveEventManager::AssertGeometry();
274799f3 116 AliFMDGeometry* geom = AliFMDGeometry::Instance();
117 geom->Init();
118 geom->InitTransformations();
119
120 AliFMDParameters* pars = AliFMDParameters::Instance();
121 // pars->UseRcuTrailer(!old);
122 pars->UseCompleteHeader(old);
123 pars->SetSampleRate(4);
124 pars->Init(kFALSE, 0);
125
126 // Get shapes
127 TGeoShape* inner = static_cast<TGeoShape*>(gGeoManager->GetListOfShapes()
128 ->FindObject("FISE"));
129 if (!inner) throw TEveException("Shape of inner type sensors not found");
130 TGeoShape* outer = static_cast<TGeoShape*>(gGeoManager->GetListOfShapes()
131 ->FindObject("FOSE"));
132 if (!outer) throw TEveException("Shape of outer type sensors not found");
133
134 // Emulate reference counting
135 inner->SetUniqueID(1000);
136 outer->SetUniqueID(1000);
137
138 // Loop over detectors
139 for (UShort_t d = 1; d <= 3; d++) {
140 AliFMDDetector* detector = geom->GetDetector(d);
141 if (!detector) continue;
142 TEveElementList* ed = new TEveElementList(Form(kDetector,
143 detector->GetId()));
144 AddElement(ed);
145 ed->IncDenyDestroy();
146 ed->SetUserData(detector);
147
148 // Loop over rings
149 Char_t rings[] = { 'I', 'O', 0 };
150 Char_t* pr = &(rings[0]);
151 while (*pr) {
152 AliFMDRing* ring = detector->GetRing(*pr);
153 pr++;
154 if (!ring) continue;
155 TEveElementList* er = new TEveElementList(Form(kRing,
156 detector->GetId(),
157 ring->GetId()));
158 ed->AddElement(er);
159 er->IncDenyDestroy();
160 er->SetUserData(ring);
161
162 // UShort_t nsec = ring->GetNSectors();
163 // UShort_t nstr = ring->GetNStrips();
164 UShort_t nmod = ring->GetNModules();
165 // Loop over modules
166 for (UShort_t m = 0; m < nmod; m++) {
167 TEveGeoShape* em = new TEveGeoShape(Form(kModule,
168 detector->GetId(),
169 ring->GetId(),
170 2*m, 2*m+1));
171 er->AddElement(em);
172 em->SetTransMatrix(*(detector->FindTransform(ring->GetId(), 2*m)));
173 em->SetShape(ring->GetId() == 'I' ? inner : outer);
fbc350a3 174 em->SetMainColor(kRed);
274799f3 175 em->SetMainTransparency(32);
176 em->IncDenyDestroy();
177
178#if 0
179 for (UShort_t s = 2*m; s < 2*m+2 && s < nsec; s++) {
180 TEveDigitSet* eb = MakeDigitSet(Form(kSector,
181 detector->GetId(),
182 ring->GetId(), s), nstr);
183 em->AddElement(eb);
184 eb->SetEmitSignals(kFALSE);
185 eb->SetPickable(kTRUE);
186 // eb->SetOwnIds(kTRUE);
187 } // for (UShort_t s ...)
188#endif
189 } // for (UShort_t m ...)
190 } // while (pr)
191 } // for (UShort_t d ...)
192}
193
194//____________________________________________________________________
195AliEveFMDLoader::~AliEveFMDLoader()
196{
197 // Destructor
198 AliWarning("AliEveFMDLoader being destroyed!");
199}
200
201//____________________________________________________________________
202Int_t
203AliEveFMDLoader::RemoveFromListTrees(TEveElement* el)
204{
205 // Called when the element should be removed from the list. We
206 // overload this to allow clearing of signals.
207 // Parameters:
208 // @param el Tree to remove from.
209
210 // Since we're most likely setting up for a new event, we clear all
211 // signals here - a little tricky, but it works(tm).
212 ClearDigitSets("All");
213
214 // Do normal TEveElement::RemoveElement
215 return TEveElementList::RemoveFromListTrees(el);
216}
217//____________________________________________________________________
218void
219AliEveFMDLoader::RemoveParent(TEveElement* el)
220{
221 // Called when the element should be removed from the list. We
222 // overload this to allow clearing of signals.
223 // Parameters:
224 // @param el Parent to remove from.
225 TEveElementList::RemoveParent(el);
226}
227
228//____________________________________________________________________
229TEveDigitSet*
230AliEveFMDLoader::MakeDigitSet(const char* name, UShort_t nstr)
231{
232 // Make a digit set. The type of digit set depends on the setting
233 // of fUseBoxDigits. If this is true, we return a TEveBoxSet,
234 // otherwise a TEveQuadSet
235 // Parameters:
236 // name The name
237 // nstr The number of strips
238 // Return
239 // newly allocated digit set
240 TEveDigitSet* ret = 0;
241 if (fUseBoxDigits) {
242 TEveBoxSet* boxes = new TEveBoxSet(name);
243 // boxes->Reset(TEveBoxSet::kBT_AABox, kFALSE, nstr);
244 boxes->Reset(TEveBoxSet::kBT_FreeBox, kFALSE, nstr);
245 ret = boxes;
246 }
247 else {
248 TEveQuadSet* quads = new TEveQuadSet(name);
249 quads->Reset(TEveQuadSet::kQT_RectangleXY, kFALSE, nstr);
250 ret = quads;
251 }
252 return ret;
253}
254
255//____________________________________________________________________
256void
257AliEveFMDLoader::ClearDigitSets(const char* type)
258{
259 // Clear signals of some type.
260 // Parameters:
261 // @param type Type of signals to clear
262 // Type can be one of
263 // - All All signals
264 // - Hits Hits
265 // - Digits Digits
266 // - Raw Raw
267 // - ESD ESD
268 TString stype(type);
269
270 for (TEveElement::List_i di = BeginChildren();
271 di != EndChildren(); ++di) {
272 for (TEveElement::List_i ri = (*di)->BeginChildren();
273 ri != (*di)->EndChildren(); ++ri) {
274 for (TEveElement::List_i mi = (*ri)->BeginChildren();
275 mi != (*ri)->EndChildren(); ++mi) {
276 if (stype == "All") {
277 (*mi)->RemoveElements();
278 continue;
279 }
280 for (TEveElement::List_i si = (*mi)->BeginChildren();
281 si != (*mi)->EndChildren(); ++si) {
282 TEveDigitSet* signals = static_cast<TEveDigitSet*>((*si));
283 if (!signals) continue;
284 TString s(signals->GetName());
285 if (!s.Contains(type)) continue;
286 (*mi)->RemoveElement(signals);
287 }
288 }
289 }
290 }
291}
292
293//____________________________________________________________________
294TEveDigitSet*
295AliEveFMDLoader::FindDigitSet(const char* t, UShort_t d, Char_t r, UShort_t s)
296{
297 // Find a digit set corresponding to the passed parameters. If it
298 // is not found, one is created
299 // Parameters:
300 // @param type Type of data
301 // @param d Detector
302 // @param r Ring
303 // @param s Sector
304 // @return a digit set
305 TEveElement* detector = FindChild(Form(kDetector, d));
306 if (!detector) {
307 AliError(Form("Detector %s not found", Form(kDetector, d)));
308 return 0;
309 }
310
311 TEveElement* ring = detector->FindChild(Form(kRing, d, r));
312 if (!ring) {
313 AliError(Form("Ring %s not found", Form(kRing, d, r)));
314 return 0;
315 }
316
317 Int_t mod = 2*(s/2);
318 TEveElement* module = ring->FindChild(Form(kModule, d, r, mod, mod+1));
319 if (!module) {
320 AliError(Form("Module %s not found", Form(kModule, d, r, s, s+1)));
321 return 0;
322 }
323
324 TEveElement* sector = module->FindChild(Form(kSector, d, r, s, t));
325 TEveDigitSet* signal = static_cast<TEveDigitSet*>(sector);
326 if (!sector) {
327 AliFMDRing* rng = AliFMDGeometry::Instance()->GetRing(r);
328 signal = MakeDigitSet(Form(kSector, d, r, s, t), rng->GetNStrips());
329 module->AddElement(signal);
330 signal->SetEmitSignals(kFALSE);
331 signal->SetPickable(kTRUE);
332 TString st(t);
333 if (t == kHits) signal->SetPalette(&fHitPalette);
334 else if (t == kDigits) signal->SetPalette(&fDigitPalette);
335 else if (t == kRaw) signal->SetPalette(&fDigitPalette);
336 else if (t == kESD) {
337 signal->SetPalette(&fMultPalette);
338 signal->SetOwnIds(kTRUE);
339 }
340 }
341 return signal;
342}
343
344//____________________________________________________________________
345void
346AliEveFMDLoader::AddSignal(const char* t,
347 UShort_t det, Char_t rng, UShort_t sec,
348 UShort_t str, Float_t signal, Float_t min,
349 Float_t max, TObject* ref)
350{
351 // Add a signal to a digit set
352 // Parameters:
353 // @param type Type of data
354 // @param det Detector
355 // @param rng Ring
356 // @param sec Sector
357 // @param str Strip
358 // @param signal Signal value
359 // @param min Minimum of this kind of signal
360 // @param max Maximum of this kind of signal
361 // @param ref Reference object
362 AliFMDGeometry* geom = AliFMDGeometry::Instance();
363 Double_t x, y, z;
364 geom->Detector2XYZ(det, rng, sec, str, x, y, z);
365 AddSignal(t, det, rng, sec, str, x, y, z, signal, min, max, ref);
366}
367
368//____________________________________________________________________
369void
370AliEveFMDLoader::AddSignal(const char* t,
371 UShort_t det, Char_t rng, UShort_t sec,
372 UShort_t str, Double_t x, Double_t y, Double_t z,
373 Float_t signal, Float_t min, Float_t max,
374 TObject* ref)
375{
376 // Add a signal to a digit set, with known (x,y,z) coordinates
377 // (this is for hits)
378 // Parameters:
379 // @param type Type of data
380 // @param det Detector
381 // @param rng Ring
382 // @param sec Sector
383 // @param str Strip
384 // @param x X coordinate
385 // @param y Y coordinate
386 // @param z Z coordinate
387 // @param signal Signal value
388 // @param min Minimum of this kind of signal
389 // @param max Maximum of this kind of signal
390 // @param ref Reference object
391 AliFMDGeometry* geom = AliFMDGeometry::Instance();
392 AliFMDRing* ring = geom->GetRing(rng);
393 if (!ring) return;
394
395 TEveDigitSet* signals = FindDigitSet(t, det, rng, sec);
396 if (!signals) {
397 AliWarning(Form("No signal (%s) found for FMD%d%c[%02d,%03d]",
398 t, det, rng, sec, str));
399 return;
400 }
401
402 Float_t scaled = TMath::Min((signal - min) / (max - min) * 10., 10.);
403 Double_t w = 2*ring->GetPitch();
404 Int_t value = int(TMath::Nint(signal));
405 AliDebug(1, Form("New signal at FMD%d%c[%02d,%03d]=%f (v=%d, s=%f)",
406 det, rng, sec, str, signal, value, scaled));
407 AddDigit(signals, x, y, z, w, scaled, value, ref);
408}
409
410//____________________________________________________________________
411void
412AliEveFMDLoader::AddDigit(TEveDigitSet* signals,
413 Double_t x, Double_t y, Double_t z,
414 Double_t w, Float_t scaled, Int_t value,
415 TObject* ref)
416{
417 // Add a digit to a digit set.
418 // Parameters:
419 // @param signals Digit set.
420 // @param x X coordinate
421 // @param y Y coordinate
422 // @param z Z coordinate
423 // @param w strip pitch
424 // @param scaled Scaled value
425 // @param value Signal value
426 // @param ref Reference object
427 if (fUseBoxDigits) {
428 TEveBoxSet* boxes = static_cast<TEveBoxSet*>(signals);
429 Float_t zc = (z > 0 ? -1 : 1) * scaled + z;
430 Float_t vs[] = { -w, -5*w, zc-scaled, // Lower back left
431 +w, -5*w, zc-scaled, // Lower back right
432 +w, +5*w, zc-scaled, // Lower front right
433 -w, +5*w, zc-scaled, // Lower front left
434 -w, -5*w, zc+scaled, // Upper back left
435 +w, -5*w, zc+scaled, // Upper back right
436 +w, +5*w, zc+scaled, // Upper front right
437 -w, +5*w, zc+scaled }; // Upper front left
438 Float_t ang = TMath::ATan2(y,x);
439 for (size_t i = 0; i < 8; i++) {
440 Float_t bx = vs[3*i+0];
441 Float_t by = vs[3*i+1];
442 Float_t ca = TMath::Cos(ang);
443 Float_t sa = TMath::Sin(ang);
444 vs[3*i+0] = bx * ca - by * sa + x;
445 vs[3*i+1] = bx * sa + by * ca + y;
446 }
447 // boxes->AddBox(x, y, (z > 0 ? -scaled : 0) + z , 5*w, w, scaled);
448 boxes->AddBox(vs);
449 boxes->DigitValue(value);
450 if (ref) boxes->DigitId(ref);
451 }
452 else {
453 TEveQuadSet* quads = static_cast<TEveQuadSet*>(signals);
454 quads->AddQuad(x,y,z,w,w);
455 quads->QuadValue(value);
456 if (ref) quads->QuadId(ref);
457 }
458}
459
460
461//____________________________________________________________________
462void
463AliEveFMDLoader::CheckAdd()
464{
465 // check if we shoul re-add ourselves to the current event node
466 TEveElement* event = gEve->GetCurrentEvent();
467 if (event && event->FindChild(GetName())) return;
468 gEve->AddElement(this);
469}
470
471
472//____________________________________________________________________
473void
474AliEveFMDLoader::LoadHits()
475{
476 // Load and display hits
477 ClearDigitSets(kHits);
478
479 AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
480 if (!rl) {
481 AliError("No run loader");
482 return;
483 }
484
485 rl->LoadHits("FMD");
486 TTree* ht = rl->GetTreeH("FMD", false);
487 if (!ht) {
488 AliError("No FMD tree");
489 return;
490 }
491
492 TClonesArray* hits = &fHitCache;
493 fHitCache.Clear();
494 ht->SetBranchAddress("FMD", &hits);
495
496 Float_t min = fHitPalette.GetMinVal();
497 Float_t max = fHitPalette.GetMaxVal();
498
499 Int_t nTracks = ht->GetEntriesFast();
500 for (Int_t i = 0; i < nTracks; i++) {
501 Int_t hitRead = ht->GetEntry(i);
502 if (hitRead <= 0) continue;
503
504 Int_t nHit = hits->GetEntriesFast();
505 if (nHit <= 0) continue;
506
507 for (Int_t j = 0; j < nHit; j++) {
508 AliFMDHit* hit = static_cast<AliFMDHit*>(hits->At(j));
509 if (!hit) continue;
510
511 AddSignal(kHits,
512 hit->Detector(), hit->Ring(), hit->Sector(), hit->Strip(),
513 hit->X(), hit->Y(), hit->Z(), int(hit->Edep()*1000),
514 min, max, hit);
515 }
516 }
517 CheckAdd();
518}
519
520//____________________________________________________________________
521void
522AliEveFMDLoader::DoLoadDigits(const char* t, TClonesArray* digits)
523{
524 // Do the actual display of digits
525 // Parameters:
526 // @param type What to show
527 // @param digits The digits
528 Float_t min = fDigitPalette.GetMinVal();
529 Float_t max = fDigitPalette.GetMaxVal();
530
531 Int_t n = digits->GetEntriesFast();
532 for (Int_t i = 0; i < n; i++) {
533 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
534 if (!digit) return;
535 AddSignal(t, digit->Detector(), digit->Ring(), digit->Sector(),
536 digit->Strip(), digit->Counts(), min, max, digit);
537 }
538 CheckAdd();
539}
540
541//____________________________________________________________________
542void
543AliEveFMDLoader::LoadDigits()
544{
545 // Load and display simulated digits
546 ClearDigitSets(kDigits);
547
548 AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
549 if (!rl) {
550 AliError("No run-loader");
551 return;
552 }
553
554 rl->LoadDigits("FMD");
555 TTree* dt = rl->GetTreeD("FMD", false);
556 if (!dt) {
557 AliError("No FMD tree");
558 return;
559 }
560
561 TClonesArray* digits = &fDigitCache;
562 fDigitCache.Clear();
563 dt->SetBranchAddress("FMD", &digits);
564
565 Int_t read = dt->GetEntry(0);
566 if (read <= 0) {
567 AliWarning("Nothing read");
568 return;
569 }
570 DoLoadDigits(kDigits, digits);
571}
572
573
574//____________________________________________________________________
575void
576AliEveFMDLoader::LoadRaw()
577{
578 // Load and display raw digits
579 ClearDigitSets(kRaw);
580
581 AliRawReader* rr = AliEveEventManager::AssertRawReader();
582 if (!rr) {
583 AliError("No raw-reader");
584 return;
585 }
586 rr->Reset();
587 std::cout<<"Now in event # " << *(rr->GetEventId()) << std::endl;
588 AliFMDRawReader* fr = new AliFMDRawReader(rr, 0);
589 TClonesArray* digits = &fRawCache;
590 fRawCache.Clear();
591
592 fr->ReadAdcs(digits);
593
594 DoLoadDigits(kRaw, digits);
595}
596
597//____________________________________________________________________
598void
599AliEveFMDLoader::LoadESD()
600{
601 // Load and display ESD information
602 ClearDigitSets(kESD);
603
604 AliESDEvent* esd = AliEveEventManager::AssertESD();
605 if (!esd) {
606 AliError("No ESD");
607 return;
608 }
609
610 AliESDFMD* fmd = esd->GetFMDData();
611 if (!fmd) {
612 AliError("No FMD ESD data");
613 return;
614 }
615
616 Float_t min = fMultPalette.GetMinVal();
617 Float_t max = fMultPalette.GetMaxVal();
618
619 for (UShort_t det = 1; det <= 3; det++) {
620 Char_t rings[] = { 'I', (det == 1 ? '\0' : 'O'), '\0' };
621 for (Char_t* rng = rings; *rng != '\0'; rng++) {
622 UShort_t nsec = (*rng == 'I' ? 20 : 40);
623 UShort_t nstr = (*rng == 'I' ? 512 : 256);
624 for (UShort_t sec = 0; sec < nsec; sec++) {
625 for (UShort_t str = 0; str < nstr; str++) {
626 Float_t mult = fmd->Multiplicity(det,*rng,sec,str);
627 if (mult == AliESDFMD::kInvalidMult) continue;
628 Float_t eta = fmd->Eta(det,*rng,sec,str);
629 AddSignal(kESD, det, *rng, sec, str, mult, min, max,
630 new TNamed(Form("FMD%d%c[%02d,%03d]", det, *rng, sec, str),
631 Form("Mch=%f, eta=%f", mult, eta)));
632 }
633 }
634 }
635 }
636 CheckAdd();
637}
638
639
640
641
642
643
644//____________________________________________________________________
645//
646// EOF
647//