AliRoot Core  v5-06-15 (45dab64)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliMUONPainterMatrix.cxx
Go to the documentation of this file.
1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15 
16 // $Id$
17 
18 #include "AliMUONPainterMatrix.h"
19 
20 #include "AliLog.h"
21 #include "AliMUONPainterGroup.h"
22 #include "AliMUONPainterHelper.h"
23 #include "AliMUONVPainter.h"
24 #include "AliMUONVTrackerData.h"
25 #include "TCanvas.h"
26 #include "TGClient.h"
27 #include "TPaveLabel.h"
28 #include <Riostream.h>
29 #include <TBox.h>
30 #include <TMath.h>
31 #include <TObjArray.h>
32 #include <TObjString.h>
33 #include <TROOT.h>
34 #include <TVirtualPad.h>
35 #include <float.h>
36 #include "AliMUONPainterEnv.h"
37 #include <cassert>
38 
44 
48 
49 //_____________________________________________________________________________
50 AliMUONPainterMatrix::AliMUONPainterMatrix(const char* name, Int_t nx, Int_t ny)
51 : TObject(),
52  fBasename(name),
53  fWhatname(""),
54  fNx(nx),
55  fNy(ny),
56  fPainters(new TObjArray(fNx*fNy)),
57  fAttributes(),
58  fName()
59 {
61 
62  fPainters->SetOwner(kTRUE);
63  if ( fNx*fNy > 1 )
64  {
65  fAttributes.SetSingle(kFALSE);
66  }
67  SetName();
68 }
69 
70 //_____________________________________________________________________________
72 {
74  delete fPainters;
75 }
76 
77 //_____________________________________________________________________________
78 void
80 {
82  fPainters->AddLast(painter);
84 }
85 
86 
87 //_____________________________________________________________________________
90 {
92 
94 
95  for ( Int_t i = 0; i < Size(); ++i )
96  {
97  AliMUONVPainter* oldPainter = Painter(i);
98 
99  AliMUONVPainter* newPainter(0x0);
100 
101  newPainter = AliMUONVPainter::CreatePainter(oldPainter->ClassName(),
102  attributes,
103  oldPainter->ID0(),
104  oldPainter->ID1());
105 
106  if (newPainter)
107  {
108  newPainter->UpdateGroupsFrom(*(oldPainter->Master()));
109  clone->Adopt(newPainter);
110  }
111  else
112  {
113  AliError(Form("Failed to create painter of class %s ID0 %d ID1 %d",
114  oldPainter->ClassName(),
115  oldPainter->ID0(),
116  oldPainter->ID1()));
117  }
118  }
119 
120  return clone;
121 }
122 
123 //_____________________________________________________________________________
124 void
126 {
128 
129  Double_t dataMin(FLT_MAX);
130  Double_t dataMax(-FLT_MAX);
131  Bool_t atLeastOnePlotter(kFALSE);
132 
133  for ( Int_t i = 0; i < Size(); ++i )
134  {
135  AliMUONVPainter* p = Painter(i);
137 
138  Double_t min(FLT_MAX);
139  Double_t max(-FLT_MAX);
140 
141  if ( g )
142  {
143  atLeastOnePlotter = kTRUE;
144  g->ComputeDataRange(min,max);
145  if ( min <= max )
146  {
147  dataMin = TMath::Min(min,dataMin);
148  dataMax = TMath::Max(max,dataMax);
149  }
150  }
151  }
152 
153  if ( dataMin > dataMax && atLeastOnePlotter )
154  {
155  AliError(Form("data min %e > max %e : setting both to 0.0",
156  dataMin,dataMax));
157  dataMin = dataMax = 0.0;
158  }
159 
160  SetDataRange(dataMin,dataMax);
161 }
162 
163 //_____________________________________________________________________________
164 void
165 AliMUONPainterMatrix::Connect(const char* sourceMethod, const char* destClassName,
166  void* destObject, const char* destMethod)
167 {
169 
170  for ( Int_t i = 0; i < Size(); ++i )
171  {
172  Painter(i)->Connect(sourceMethod,destClassName,destObject,destMethod);
173  }
174 }
175 
176 //_____________________________________________________________________________
177 TCanvas*
178 AliMUONPainterMatrix::CreateCanvas(Int_t x, Int_t y, Int_t w, Int_t h)
179 {
198 
199  Int_t mw = ( w <= 0 ? TMath::Nint(gClient->GetDisplayWidth()*0.9) : w );
200  Int_t mh = ( h <= 0 ? TMath::Nint(gClient->GetDisplayHeight()*0.9) : h );
201 
202  TString name(GetName());
203 
204  TCanvas* d = new TCanvas(name.Data(),name.Data(),x,y,mw,mh);
205 
206  TVirtualPad* pTitle = new TPad(Form("%s-title",name.Data()),Form("%s-title",name.Data()),0,0.9,1.0,0.99);
207 
208  pTitle->Draw();
209 
210  pTitle->cd();
211 
212  TPaveLabel* text = new TPaveLabel(0,0,1,1,"");
213  text->SetFillStyle(0);
214  text->SetFillColor(0);
215  text->SetTextColor(4);
216  text->SetBorderSize(0);
217 
218  text->SetLabel(name.Data());
219 
220  text->Draw();
221 
222  d->cd();
223 
224  TVirtualPad* pMatrix = new TPad(Form("%s-matrix",name.Data()),Form("%s-matrix",name.Data()),0,0,0.9,0.89);
225 
226  pMatrix->Draw();
227  pMatrix->cd();
228 
229  Draw();
230 
231  d->cd();
232 
233  TVirtualPad* pColor = new TPad(Form("%s-color",name.Data()),Form("%s-color",name.Data()),0.91,0.01,0.99,0.89);
234 
235  pColor->Range(0,0,1,1);
236 
237  pColor->Draw();
238 
239  pColor->cd();
240 
241  Int_t ndivisions(20);
242 
243  Double_t rangeXmin(0.1);
244  Double_t rangeXmax(0.9);
245 
246  Double_t ymin, ymax;
247 
248  GetDataRange(ymin,ymax);
249 
250  Double_t min(0.0);
251  Double_t max(1.0);
252 
253  Double_t step = (max-min)/ndivisions;
254 
255  Double_t hsize = 1.0/(ndivisions+2);
256 
257  Double_t ypos = 1.0;
258 
259  for ( Int_t i = -1; i < ndivisions+1; ++i )
260  {
261  Double_t value = max - (min + step*i);
262 
263  Int_t color = AliMUONPainterHelper::Instance()->ColorFromValue(value,min,max);
264 
265  Bool_t limit(kFALSE);
266 
267  TString label;
268  TString sign;
269 
270  Double_t yvalue(0.0);
271 
272  if ( i == -1 )
273  {
274  yvalue = ymax;
275  limit = kTRUE;
276  sign = ">";
277  }
278  else if ( i == ndivisions )
279  {
280  yvalue = ymin;
281  limit = kTRUE;
282  sign = "<=";
283  }
284 
285  if (limit)
286  {
287  if ( TMath::Abs(yvalue) < 1E5 )
288  {
289  label = Form("%s %7.2f",sign.Data(),yvalue);
290  }
291  else
292  {
293  label = Form("%s %e",sign.Data(),yvalue);
294  }
295  }
296 
297  TPaveLabel* box = new TPaveLabel(rangeXmin,TMath::Max(0.001,ypos-hsize),rangeXmax,ypos,label.Data(),"");
298 
299  ypos -= hsize;
300 
301  box->SetFillColor(color);
302  box->SetTextColor( i == -1 ? 0 : 1 );
303  box->SetBorderSize(1);
304  box->SetLineColor(1);
305  box->Draw();
306  }
307 
308  d->SetEditable(kFALSE);
309 
310  return d;
311 }
312 
313 //_____________________________________________________________________________
314 void
315 AliMUONPainterMatrix::GetDataRange(Double_t& dataMin, Double_t& dataMax) const
316 {
318 
319  dataMin=FLT_MAX;
320  dataMax=-FLT_MAX;
321 
322  for ( Int_t i = 0; i < Size(); ++i )
323  {
324  AliMUONVPainter* p = Painter(i);
325  if ( p )
326  {
328  if ( g )
329  {
330  dataMin = TMath::Min(dataMin,g->DataMin());
331  dataMax = TMath::Max(dataMax,g->DataMax());
332  }
333  }
334  }
335 }
336 
337 //_____________________________________________________________________________
338 void
340 {
342 
343  types.SetOwner(kTRUE);
344  types.Clear();
345 
346  for ( Int_t i = 0; i < Size(); ++i )
347  {
348  AliMUONVPainter* p = Painter(i);
349  TObjArray ptypes;
350  p->GetTypes(ptypes);
351  TIter next(&ptypes);
352  TObject* o;
353  while ( ( o = next() ) )
354  {
355  if ( ! types.FindObject(o) )
356  {
357  types.AddLast(o->Clone());
358  }
359  }
360  }
361 }
362 
363 
364 //_____________________________________________________________________________
367 {
370  return ( group ? group->Data() : 0x0 );
371 }
372 
373 //_____________________________________________________________________________
374 TString
376 {
379  return ( group ? group->Type() : "" );
380 }
381 
382 //_____________________________________________________________________________
383 Int_t
385 {
388  return ( group ? group->DataIndex() : -1 );
389 }
390 
391 //_____________________________________________________________________________
392 void
394 {
396 
397  if (!gPad)
398  {
399  gROOT->MakeDefCanvas();
400  }
401 
402  TVirtualPad* pad = gPad;
403 
404  gPad->Divide(Nx(),Ny());
405 
406  for ( Int_t i = 0; i < Size(); ++i )
407  {
408  AliMUONVPainter* painter = Painter(i);
409  pad->cd(i+1);
410  painter->Draw("R");
411  }
412 
413  AppendPad("");
414 }
415 
416 //_____________________________________________________________________________
417 std::string
418 AliMUONPainterMatrix::NameIt(const char* whatname, const char* basename, const AliMUONAttPainter& att)
419 {
421  if ( strlen(whatname) > 0 )
422  {
423  return Form("%s-%s-%s",whatname,basename,att.GetName());
424  }
425  else
426  {
427  return Form("nodata-%s-%s",basename,att.GetName());
428  }
429 }
430 
431 //_____________________________________________________________________________
434 {
436 
437  if ( index <= fPainters->GetLast() )
438  {
439  return static_cast<AliMUONVPainter*>(fPainters->At(index));
440  }
441  return 0x0;
442 }
443 
444 //_____________________________________________________________________________
445 void
447 {
449  std::cout << "Whatname=" << fWhatname.Data() << " Basename=" << fBasename.Data()
450  << " Nx=" << fNx << " Ny=" << fNy << " Att=" << fAttributes.GetName() << std::endl;
451 }
452 
453 //_____________________________________________________________________________
454 void
456  Int_t indexInData)
457 {
459 
460  for ( Int_t i = 0; i < Size(); ++i )
461  {
462  AliMUONVPainter* painter = Painter(i);
463  painter->SetData(pattern,d,indexInData);
464  }
465 
466  if ( d )
467  {
468  fWhatname = Form("%s-%s",d->GetName(),d->DimensionName(indexInData).Data());
469  }
470  else
471  {
472  fWhatname = "";
473  }
474 
475  SetName();
476 
477  // get the data source range, if any is given
478 
479  if ( d )
480  {
482 
483  TString desc = env->DataSourceDescriptor(d->GetName());
484  TString dimensionName = d->DimensionName(DataIndex());
485 
486  Double_t xmin,xmax;
487 
488  Bool_t ok = env->Ranges2DimensionRange(env->Descriptor2Ranges(desc),dimensionName.Data(),xmin,xmax);
489 
490  if (ok)
491  {
492  SetDataRange(xmin,xmax);
493  env->Save();
494  }
495  }
496 }
497 
498 //_____________________________________________________________________________
499 void
500 AliMUONPainterMatrix::SetDataRange(Double_t dataMin, Double_t dataMax)
501 {
503 
504  for ( Int_t i = 0; i < Size(); ++i )
505  {
506  AliMUONVPainter* p = Painter(i);
508  if ( g )
509  {
510  g->SetDataRange(dataMin,dataMax);
511  }
512  }
513 
514  AliMUONVTrackerData* data = Data();
515 
516  if ( data )
517  {
519 
520  TString dimensionName = data->DimensionName(DataIndex());
521  env->SetDimensionRange(data->GetName(),dimensionName,dataMin,dataMax);
522  }
523 
524 }
525 
526 //_____________________________________________________________________________
528 {
530 // fName = NameIt(fWhatname.Data(),fBasename.Data(),fAttributes);
531  fName = "nodata";
532 
533  if ( fWhatname.Length() > 0 )
534  {
535  fName = fWhatname;
536  }
537  fName += "-";
538  fName += fBasename;
539  fName += "-";
541 }
542 
543 
544 //_____________________________________________________________________________
545 //void
546 //AliMUONPainterMatrix::ChangeAttributes(const AliMUONAttPainter& attributes)
547 //{
548 // /// Change painters' attributes
549 //
550 // AliWarning("Implement me !");
551 //
552 // // for ( Int_t i = 0; i < Size(); ++i )
553 // // {
554 // // Painter(i)->SetAttributes(attributes);
555 // // }
556 //}
557 
558 
559 //_____________________________________________________________________________
560 void
561 AliMUONPainterMatrix::SetOutlined(const char* pattern, Bool_t value)
562 {
564 
565  for ( Int_t i = 0; i < Size(); ++i )
566  {
567  Painter(i)->SetOutlined(pattern,value);
568  }
569 }
570 
571 //_____________________________________________________________________________
572 void
574 {
576  for ( Int_t i = 0; i < Size(); ++i )
577  {
578  Painter(i)->SetResponder(pattern);
579  }
580 }
581 
582 //_____________________________________________________________________________
583 Int_t
585 {
587  return fPainters->GetLast()+1;
588 }
589 
590 //_____________________________________________________________________________
591 void
593 {
595 
596  Bool_t cathode0(kFALSE);
597  Bool_t cathode1(kFALSE);
598  Bool_t bending(kFALSE);
599  Bool_t nonbending(kFALSE);
600  Bool_t front(kFALSE);
601  Bool_t back(kFALSE);
602  Bool_t cathplaneexclusive(kFALSE);
603  Bool_t cathplanedisabled(kFALSE);
604 
605  for ( Int_t i = 0; i < Size(); ++i )
606  {
608 
609  if ( att.IsCathodeDefined() )
610  {
611  if ( att.IsCathode0() ) cathode0 = kTRUE;
612  if ( att.IsCathode1() ) cathode1 = kTRUE;
613  }
614 
615  if ( att.IsPlaneDefined() )
616  {
617  if ( att.IsBendingPlane() ) bending = kTRUE;
618  if ( att.IsNonBendingPlane() ) nonbending = kTRUE;
619  }
620 
621  if ( att.IsFrontView() ) front = kTRUE;
622  if ( att.IsBackView() ) back = kTRUE;
623 
624  if ( att.IsCathodeAndPlaneMutuallyExclusive() ) cathplaneexclusive = kTRUE;
625 
626  if ( att.IsCathodeAndPlaneDisabled() ) cathplanedisabled = kTRUE;
627  }
628 
629  fAttributes.SetCathode(cathode0,cathode1);
630  fAttributes.SetPlane(bending,nonbending);
631  fAttributes.SetViewPoint(front,back);
633  fAttributes.SetCathodeAndPlaneDisabled(cathplanedisabled);
634 
635  SetName();
636 }
637 
638 //_____________________________________________________________________________
641 {
643 
645 
646  for ( Int_t i = 0; i < Size() && a.IsValid(); ++i )
647  {
648  a = Painter(i)->Validate(att);
649  }
650  return a;
651 }
Int_t ColorFromValue(Double_t value, Double_t min, Double_t max) const
static AliMUONPainterHelper * Instance()
TObjArray * fPainters
painters in that matrix
Bool_t IsNonBendingPlane() const
Whether we are representing non bending plane.
virtual const char * GetName() const
Get matrix name.
Int_t Nx() const
Number of painters to arrange in x-direction.
void SetOutlined(const char *pattern, Bool_t flag)
#define TObjArray
Int_t DataIndex() const
Return the index within the data.
virtual void Draw(Option_t *opt="")
Bool_t IsPlaneDefined() const
Whether we are defined by plane.
void Draw(Option_t *opt="")
void SetResponder(const char *pattern)
virtual void SetResponder(const char *pattern)
TROOT * gROOT
Bool_t IsCathodeDefined() const
Whether we are defined by cathode.
Bool_t IsBendingPlane() const
Whether we are representing bending plane.
static std::string NameIt(const char *what, const char *basename, const AliMUONAttPainter &att)
void SetViewPoint(Bool_t front, Bool_t back)
void SetCathode(Bool_t cath0, Bool_t cath1)
void GetTypes(TObjArray &types) const
Bool_t IsCathodeAndPlaneDisabled() const
Whether cathode & plane are disabled.
virtual TString DimensionName(Int_t dim) const =0
Get the name of a given (internal) dimension.
Int_t ID0() const
Get our first ID.
void SetDataRange(Double_t min, Double_t max)
Set the data range.
static AliMUONVPainter * CreatePainter(const char *className, const AliMUONAttPainter &att, Int_t id1, Int_t id2)
Bool_t IsBackView() const
Whether the painter is to be represented from back (as seen from IP)
TString fBasename
base name of that matrix
Resource file handling.
ClassImp(TPCGenInfo)
Definition: AliTPCCmpNG.C:254
AliMUONVTrackerData * Data() const
Return the data we are plotting.
void Connect(const char *sourceMethod, const char *destClassName, void *destObject, const char *destMethod)
void ComputeDataRange()
Compute the data range for this matrix.
static TString Descriptor2Ranges(const char *dataSourceDescriptor)
TCanvas * CreateCanvas(Int_t x=0, Int_t y=0, Int_t w=-1, Int_t h=-1)
const AliMUONAttPainter & Attributes() const
Get our attributes.
const char * Type() const
Our type.
virtual AliMUONAttPainter Validate(const AliMUONAttPainter &attributes) const
Convert attributes so they are valid ones for us.
Int_t ID1() const
Get our second ID.
Base class for a graphical object representing some part of the MUON tracking system.
Bool_t IsCathodeAndPlaneMutuallyExclusive() const
Whether we can select both cathode and plane.
void SetDimensionRange(const char *dataSourceName, const char *dimensionName, Double_t xmin, Double_t xmax)
const char * Basename() const
Base name (short name)
AliMUONPainterMatrix * Clone(const AliMUONAttPainter &attributes) const
Int_t fNx
number of rows
void SetCathodeAndPlaneDisabled(Bool_t value)
void ComputeDataRange(Double_t &dataMin, Double_t &dataMax)
void SetData(const char *pattern, AliMUONVTrackerData *d, Int_t indexInData)
Bool_t IsCathode1() const
Whether we are representing cathode 1.
virtual const char * GetName() const
Return our name.
TString DataSourceDescriptor(const char *dataSourceName) const
AliMUONVPainter * Master() const
TString fWhatname
data name
void GetTypes(TObjArray &types) const
void GetDataRange(Double_t &dataMin, Double_t &dataMax) const
Get the data range for this matrix.
AliMUONPainterMatrix(const char *basename="", Int_t nx=1, Int_t ny=1)
AliMUONAttPainter Validate(const AliMUONAttPainter &att) const
Normalize attributes.
static Bool_t Ranges2DimensionRange(const char *ranges, const char *dimensionName, Double_t &xmin, Double_t &xmax)
AliMUONVTrackerData * Data() const
void Adopt(AliMUONVPainter *painter)
Adopt a painter in this matrix.
Bool_t IsFrontView() const
Whether the painter is to be represented from front (as seen from IP)
Bool_t IsValid() const
Whether we are valid.
Bool_t IsCathode0() const
Whether we are representing cathode 0.
AliMUONPainterEnv * Env()
Return the environment.
void SetDataRange(Double_t min, Double_t max)
Force a given data range for all painter groups belonging to this matrix.
void UpdateGroupsFrom(const AliMUONVPainter &painter)
Basic attributes shared by all painters.
Double_t DataMin() const
Min data we are plotting.
void SetCathodeAndPlaneMutuallyExclusive(Bool_t value)
void SetData(const char *pattern, AliMUONVTrackerData *data, Int_t dataIndex)
A group of AliMUONVPainter.
void SetPlane(Bool_t bending, Bool_t nonBending)
AliMUONPainterGroup * PlotterGroup() const
Return the plotter group.
Int_t Ny() const
Number of painters to arrange in y-direction.
Base class for MUON data that can be presented at different levels in the hierarchy of the MUON syste...
AliMUONVPainter * Painter(Int_t index) const
Get a painter.
A matrix of AliMUONVPainter.
void SetOutlined(const char *pattern, Bool_t value)
Double_t DataMax() const
Max data we are plotting.
Int_t Size() const
Number of painters (should be <= Nx*Ny)
Int_t fNy
number of columns
void Print(Option_t *opt="") const
Printout.
AliMUONAttPainter fAttributes
attributes of our painter(s)