Logo ROOT   6.08/07
Reference Guide
TMarker.cxx
Go to the documentation of this file.
1 // @(#)root/graf:$Id$
2 // Author: Rene Brun 12/05/95
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 #include <stdlib.h>
13 
14 #include "Riostream.h"
15 #include "TROOT.h"
16 #include "TVirtualPad.h"
17 #include "TMarker.h"
18 #include "TVirtualX.h"
19 #include "TMath.h"
20 #include "TPoint.h"
21 #include "TText.h"
22 #include "TClass.h"
23 #include "TPoint.h"
24 
26 
27 
28 /** \class TMarker
29 \ingroup BasicGraphics
30 
31 Manages Markers. Marker attributes are managed by TAttMarker.
32 */
33 
34 ////////////////////////////////////////////////////////////////////////////////
35 /// Marker default constructor.
36 
38 {
39  fX = 0;
40  fY = 0;
41 }
42 
43 ////////////////////////////////////////////////////////////////////////////////
44 /// Marker normal constructor.
45 
47  :TObject(), TAttMarker()
48 {
49  fX = x;
50  fY = y;
51  fMarkerStyle = marker;
52 }
53 
54 ////////////////////////////////////////////////////////////////////////////////
55 /// Marker default destructor.
56 
58 {
59 }
60 
61 ////////////////////////////////////////////////////////////////////////////////
62 /// Marker copy constructor.
63 
64 TMarker::TMarker(const TMarker &marker) : TObject(marker), TAttMarker(marker), TAttBBox2D(marker)
65 {
66  fX = 0;
67  fY = 0;
68  ((TMarker&)marker).Copy(*this);
69 }
70 
71 ////////////////////////////////////////////////////////////////////////////////
72 /// Copy this marker to marker.
73 
74 void TMarker::Copy(TObject &obj) const
75 {
76  TObject::Copy(obj);
77  TAttMarker::Copy(((TMarker&)obj));
78  ((TMarker&)obj).fX = fX;
79  ((TMarker&)obj).fY = fY;
80 }
81 
82 ////////////////////////////////////////////////////////////////////////////////
83 /// Display the table of markers with their numbers.
84 
86 {
87  TMarker *marker = new TMarker();
88  marker->SetMarkerSize(3);
89  TText *text = new TText();
90  text->SetTextFont(62);
91  text->SetTextAlign(22);
92  text->SetTextSize(0.1);
93  char atext[] = " ";
94  Double_t x = 0;
95  Double_t dx = 1/16.0;
96  for (Int_t i=1;i<16;i++) {
97  x += dx;
98  snprintf(atext,7,"%d",i);
99  marker->SetMarkerStyle(i);
100  marker->DrawMarker(x,.35);
101  text->DrawText(x,.17,atext);
102  snprintf(atext,7,"%d",i+19);
103  marker->SetMarkerStyle(i+19);
104  marker->DrawMarker(x,.8);
105  text->DrawText(x,.62,atext);
106  }
107  delete marker;
108  delete text;
109 }
110 
111 ////////////////////////////////////////////////////////////////////////////////
112 /// Compute distance from point px,py to a marker.
113 ///
114 /// Compute the closest distance of approach from point px,py to this marker.
115 /// The distance is computed in pixels units.
116 
118 {
119  Int_t pxm, pym;
120  if (TestBit(kMarkerNDC)) {
121  pxm = gPad->UtoPixel(fX);
122  pym = gPad->VtoPixel(fY);
123  } else {
124  pxm = gPad->XtoAbsPixel(gPad->XtoPad(fX));
125  pym = gPad->YtoAbsPixel(gPad->YtoPad(fY));
126  }
127  Int_t dist = (Int_t)TMath::Sqrt((px-pxm)*(px-pxm) + (py-pym)*(py-pym));
128 
129  //marker size = 1 is about 8 pixels
130  Int_t markerRadius = Int_t(4*fMarkerSize);
131  if (dist <= markerRadius) return 0;
132  if (dist > markerRadius+3) return 999;
133  return dist;
134 }
135 
136 ////////////////////////////////////////////////////////////////////////////////
137 /// Draw this marker with its current attributes.
138 
139 void TMarker::Draw(Option_t *option)
140 {
141  AppendPad(option);
142 
143 }
144 
145 ////////////////////////////////////////////////////////////////////////////////
146 /// Draw this marker with new coordinates.
147 
149 {
150  TMarker *newmarker = new TMarker(x, y, 1);
151  TAttMarker::Copy(*newmarker);
152  newmarker->SetBit(kCanDelete);
153  newmarker->AppendPad();
154 }
155 
156 ////////////////////////////////////////////////////////////////////////////////
157 /// Execute action corresponding to one event.
158 ///
159 /// This member function is called when a marker is clicked with the locator
160 ///
161 /// If Left button is clicked on a marker, the marker is moved to
162 /// a new position when the mouse button is released.
163 
165 {
166  if (!gPad) return;
167 
168  TPoint p;
169  static Int_t pxold, pyold;
170  static Bool_t ndcsav;
171  Double_t dpx, dpy, xp1,yp1;
172  Bool_t opaque = gPad->OpaqueMoving();
173 
174  if (!gPad->IsEditable()) return;
175 
176  switch (event) {
177 
178  case kButton1Down:
179  ndcsav = TestBit(kMarkerNDC);
180  if (!opaque) {
181  gVirtualX->SetTextColor(-1); // invalidate current text color (use xor mode)
182  TAttMarker::Modify(); //Change marker attributes only if necessary
183  }
184  // No break !!!
185 
186  case kMouseMotion:
187  pxold = px; pyold = py;
188  gPad->SetCursor(kMove);
189  break;
190 
191  case kButton1Motion:
192  p.fX = pxold; p.fY = pyold;
193  if (!opaque) gVirtualX->DrawPolyMarker(1, &p);
194  p.fX = px; p.fY = py;
195  if (!opaque) gVirtualX->DrawPolyMarker(1, &p);
196  pxold = px; pyold = py;
197  if (opaque) {
198  if (ndcsav) this->SetNDC(kFALSE);
199  this->SetX(gPad->PadtoX(gPad->AbsPixeltoX(px)));
200  this->SetY(gPad->PadtoY(gPad->AbsPixeltoY(py)));
201  gPad->ShowGuidelines(this, event, 'i', true);
202  gPad->Modified(kTRUE);
203  gPad->Update();
204  }
205  break;
206 
207  case kButton1Up:
208  if (opaque) {
209  if (ndcsav && !this->TestBit(kMarkerNDC)) {
210  this->SetX((fX - gPad->GetX1())/(gPad->GetX2()-gPad->GetX1()));
211  this->SetY((fY - gPad->GetY1())/(gPad->GetY2()-gPad->GetY1()));
212  this->SetNDC();
213  }
214  gPad->ShowGuidelines(this, event);
215  } else {
216  if (TestBit(kMarkerNDC)) {
217  dpx = gPad->GetX2() - gPad->GetX1();
218  dpy = gPad->GetY2() - gPad->GetY1();
219  xp1 = gPad->GetX1();
220  yp1 = gPad->GetY1();
221  fX = (gPad->AbsPixeltoX(pxold)-xp1)/dpx;
222  fY = (gPad->AbsPixeltoY(pyold)-yp1)/dpy;
223  } else {
224  fX = gPad->PadtoX(gPad->AbsPixeltoX(px));
225  fY = gPad->PadtoY(gPad->AbsPixeltoY(py));
226  }
227  gPad->Modified(kTRUE);
228  gPad->Update();
229  gVirtualX->SetTextColor(-1);
230  }
231  break;
232  }
233 }
234 
235 ////////////////////////////////////////////////////////////////////////////////
236 /// List this marker with its attributes.
237 
238 void TMarker::ls(Option_t *) const
239 {
241  printf("Marker X=%f Y=%f marker type=%d\n",fX,fY,fMarkerStyle);
242 }
243 
244 ////////////////////////////////////////////////////////////////////////////////
245 /// Paint this marker with its current attributes.
246 
248 {
249  if (TestBit(kMarkerNDC)) {
250  Double_t u = gPad->GetX1() + fX*(gPad->GetX2()-gPad->GetX1());
251  Double_t v = gPad->GetY1() + fY*(gPad->GetY2()-gPad->GetY1());
252  PaintMarker(u,v);
253  } else {
254  PaintMarker(gPad->XtoPad(fX),gPad->YtoPad(fY));
255  }
256 }
257 
258 ////////////////////////////////////////////////////////////////////////////////
259 /// Draw this marker with new coordinates.
260 
262 {
263  TAttMarker::Modify(); //Change line attributes only if necessary
264  gPad->PaintPolyMarker(-1,&x,&y,"");
265 }
266 
267 ////////////////////////////////////////////////////////////////////////////////
268 /// Draw this marker with new coordinates in NDC.
269 
271 {
272 }
273 
274 ////////////////////////////////////////////////////////////////////////////////
275 /// Dump this marker with its attributes.
276 
278 {
279  printf("Marker X=%f Y=%f",fX,fY);
280  if (GetMarkerColor() != 1) printf(" Color=%d",GetMarkerColor());
281  if (GetMarkerStyle() != 1) printf(" MarkerStyle=%d",GetMarkerStyle());
282  if (GetMarkerSize() != 1) printf(" MarkerSize=%f",GetMarkerSize());
283  printf("\n");
284 }
285 
286 ////////////////////////////////////////////////////////////////////////////////
287 /// Save primitive as a C++ statement(s) on output stream out
288 
289 void TMarker::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
290 {
291  if (gROOT->ClassSaved(TMarker::Class())) {
292  out<<" ";
293  } else {
294  out<<" TMarker *";
295  }
296  out<<"marker = new TMarker("<<fX<<","<<fY<<","<<fMarkerStyle<<");"<<std::endl;
297 
298  SaveMarkerAttributes(out,"marker",1,1,1);
299 
300  out<<" marker->Draw();"<<std::endl;
301 }
302 
303 ////////////////////////////////////////////////////////////////////////////////
304 /// Set NDC mode on if isNDC = kTRUE, off otherwise
305 
307 {
309  if (isNDC) SetBit(kMarkerNDC);
310 }
311 
312 ////////////////////////////////////////////////////////////////////////////////
313 /// Stream an object of class TMarker.
314 
315 void TMarker::Streamer(TBuffer &R__b)
316 {
317  if (R__b.IsReading()) {
318  UInt_t R__s, R__c;
319  Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
320  if (R__v > 1) {
321  R__b.ReadClassBuffer(TMarker::Class(), this, R__v, R__s, R__c);
322  return;
323  }
324  //====process old versions before automatic schema evolution
325  TObject::Streamer(R__b);
326  TAttMarker::Streamer(R__b);
327  Float_t x,y;
328  R__b >> x; fX = x;
329  R__b >> y; fY = y;
330  //====end of old versions
331 
332  } else {
333  R__b.WriteClassBuffer(TMarker::Class(),this);
334  }
335 }
336 
337 ////////////////////////////////////////////////////////////////////////////////
338 /// Return the bounding Box of the Line
339 
341 {
342  Double_t size = this->GetMarkerSize();
343 
344  Rectangle_t BBox;
345  BBox.fX = gPad->XtoPixel(fX)+(Int_t)(2*size);
346  BBox.fY = gPad->YtoPixel(fY)-(Int_t)(2*size);
347  BBox.fWidth = 2*size;
348  BBox.fHeight = 2*size;
349  return (BBox);
350 }
351 
352 ////////////////////////////////////////////////////////////////////////////////
353 /// Return the center of the BoundingBox as TPoint in pixels
354 
356 {
357  TPoint p;
358  p.SetX(gPad->XtoPixel(fX));
359  p.SetY(gPad->YtoPixel(fY));
360  return(p);
361 }
362 
363 ////////////////////////////////////////////////////////////////////////////////
364 /// Set center of the BoundingBox
365 
367 {
368  fX = gPad->PixeltoX(p.GetX());
369  fY = gPad->PixeltoY(p.GetY() - gPad->VtoPixel(0));
370 }
371 
372 ////////////////////////////////////////////////////////////////////////////////
373 /// Set X coordinate of the center of the BoundingBox
374 
376 {
377  fX = gPad->PixeltoX(x);
378 }
379 
380 ////////////////////////////////////////////////////////////////////////////////
381 /// Set Y coordinate of the center of the BoundingBox
382 
384 {
385  fY = gPad->PixeltoY(y - gPad->VtoPixel(0));
386 }
387 
388 ////////////////////////////////////////////////////////////////////////////////
389 /// Set left hand side of BoundingBox to a value
390 /// (resize in x direction on left)
391 
393 {
394  Double_t size = this->GetMarkerSize();
395  fX = gPad->PixeltoX(x + (Int_t)size);
396 }
397 
398 ////////////////////////////////////////////////////////////////////////////////
399 /// Set right hand side of BoundingBox to a value
400 /// (resize in x direction on right)
401 
403 {
404  Double_t size = this->GetMarkerSize();
405  fX = gPad->PixeltoX(x - (Int_t)size);
406 }
407 
408 ////////////////////////////////////////////////////////////////////////////////
409 /// Set top of BoundingBox to a value (resize in y direction on top)
410 
412 {
413  Double_t size = this->GetMarkerSize();
414  fY = gPad->PixeltoY(y - (Int_t)size - gPad->VtoPixel(0));
415 }
416 
417 ////////////////////////////////////////////////////////////////////////////////
418 /// Set bottom of BoundingBox to a value
419 /// (resize in y direction on bottom)
420 
422 {
423  Double_t size = this->GetMarkerSize();
424  fY = gPad->PixeltoY(y + (Int_t)size - gPad->VtoPixel(0));
425 }
UShort_t fWidth
Definition: GuiTypes.h:364
Bool_t IsReading() const
Definition: TBuffer.h:83
virtual void DrawMarker(Double_t x, Double_t y)
Draw this marker with new coordinates.
Definition: TMarker.cxx:148
double dist(Rotation3D const &r1, Rotation3D const &r2)
Definition: 3DDistances.cxx:48
virtual void SetX(Double_t x)
Definition: TMarker.h:62
void SetX(SCoord_t x)
Definition: TPoint.h:51
Short_t fY
Definition: GuiTypes.h:363
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
TMarker()
Marker default constructor.
Definition: TMarker.cxx:37
UShort_t fHeight
Definition: GuiTypes.h:364
virtual void SetNDC(Bool_t isNDC=kTRUE)
Set NDC mode on if isNDC = kTRUE, off otherwise.
Definition: TMarker.cxx:306
virtual void Print(Option_t *option="") const
Dump this marker with its attributes.
Definition: TMarker.cxx:277
short Version_t
Definition: RtypesCore.h:61
virtual void SetBBoxCenter(const TPoint &p)
Set center of the BoundingBox.
Definition: TMarker.cxx:366
float Float_t
Definition: RtypesCore.h:53
const char Option_t
Definition: RtypesCore.h:62
virtual void SetBBoxCenterY(const Int_t y)
Set Y coordinate of the center of the BoundingBox.
Definition: TMarker.cxx:383
Double_t fX
X position of marker (left,center,etc..)
Definition: TMarker.h:34
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:157
virtual void Draw(Option_t *option="")
Draw this marker with its current attributes.
Definition: TMarker.cxx:139
SCoord_t fX
Definition: TPoint.h:37
Buffer base class used for serializing objects.
Definition: TBuffer.h:42
Size_t fMarkerSize
Marker size.
Definition: TAttMarker.h:29
virtual void SetBBoxX1(const Int_t x)
Set left hand side of BoundingBox to a value (resize in x direction on left)
Definition: TMarker.cxx:392
#define gROOT
Definition: TROOT.h:364
SCoord_t fY
Definition: TPoint.h:38
Manages Markers.
Definition: TMarker.h:31
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
SCoord_t GetY() const
Definition: TPoint.h:50
void Copy(TAttMarker &attmarker) const
Copy this marker attributes to a new TAttMarker.
Definition: TAttMarker.cxx:194
SCoord_t GetX() const
Definition: TPoint.h:49
const char * Class
Definition: TXMLSetup.cxx:64
void SetY(SCoord_t y)
Definition: TPoint.h:52
virtual void SetBBoxCenterX(const Int_t x)
Set X coordinate of the center of the BoundingBox.
Definition: TMarker.cxx:375
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:739
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive as a C++ statement(s) on output stream out.
Definition: TMarker.cxx:289
if object in a list can be deleted
Definition: TObject.h:56
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:165
virtual Style_t GetMarkerStyle() const
Return the marker style.
Definition: TAttMarker.h:37
Marker Attributes class.
Definition: TAttMarker.h:24
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition: TAttText.h:51
Double_t x[n]
Definition: legend1.C:17
virtual void Copy(TObject &object) const
Copy this to obj.
Definition: TObject.cxx:123
static void DisplayMarkerTypes()
Display the table of markers with their numbers.
Definition: TMarker.cxx:85
virtual Size_t GetMarkerSize() const
Return the marker size.
Definition: TAttMarker.h:38
Base class for several text objects.
Definition: TText.h:33
virtual void SaveMarkerAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t sizdef=1)
Save line attributes as C++ statement(s) on output stream out.
Definition: TAttMarker.cxx:229
virtual void PaintMarker(Double_t x, Double_t y)
Draw this marker with new coordinates.
Definition: TMarker.cxx:261
Short_t fX
Definition: GuiTypes.h:363
virtual void SetTextAlign(Short_t align=11)
Set the text alignment.
Definition: TAttText.h:47
Definition: TPoint.h:33
Style_t fMarkerStyle
Marker style.
Definition: TAttMarker.h:28
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
Definition: TMarker.cxx:164
virtual void PaintMarkerNDC(Double_t u, Double_t v)
Draw this marker with new coordinates in NDC.
Definition: TMarker.cxx:270
virtual void Paint(Option_t *option="")
Paint this marker with its current attributes.
Definition: TMarker.cxx:247
virtual void Modify()
Change current marker attributes if necessary.
Definition: TAttMarker.cxx:204
SVector< double, 2 > v
Definition: Dict.h:5
virtual void SetBBoxX2(const Int_t x)
Set right hand side of BoundingBox to a value (resize in x direction on right)
Definition: TMarker.cxx:402
virtual Rectangle_t GetBBox()
Return the bounding Box of the Line.
Definition: TMarker.cxx:340
unsigned int UInt_t
Definition: RtypesCore.h:42
virtual void SetBBoxY1(const Int_t y)
Set top of BoundingBox to a value (resize in y direction on top)
Definition: TMarker.cxx:411
static void IndentLevel()
Functions used by ls() to indent an object hierarchy.
Definition: TROOT.cxx:2618
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition: TAttMarker.h:45
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition: TAttMarker.h:46
#define gVirtualX
Definition: TVirtualX.h:362
virtual ~TMarker()
Marker default destructor.
Definition: TMarker.cxx:57
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
#define ClassImp(name)
Definition: Rtypes.h:279
double Double_t
Definition: RtypesCore.h:55
void Copy(TObject &marker) const
Copy this marker to marker.
Definition: TMarker.cxx:74
TText * text
Double_t y[n]
Definition: legend1.C:17
virtual void SetY(Double_t y)
Definition: TMarker.h:63
Marker position is in NDC.
Definition: TMarker.h:40
Mother of all ROOT objects.
Definition: TObject.h:37
virtual void ls(Option_t *option="") const
List this marker with its attributes.
Definition: TMarker.cxx:238
Abstract base class for elements drawn in the editor.
Definition: TAttBBox2D.h:23
virtual void SetBBoxY2(const Int_t y)
Set bottom of BoundingBox to a value (resize in y direction on bottom)
Definition: TMarker.cxx:421
#define snprintf
Definition: civetweb.c:822
virtual TText * DrawText(Double_t x, Double_t y, const char *text)
Draw this text with new coordinates.
Definition: TText.cxx:175
#define gPad
Definition: TVirtualPad.h:289
Double_t fY
Y position of marker (left,center,etc..)
Definition: TMarker.h:35
void ResetBit(UInt_t f)
Definition: TObject.h:156
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a marker.
Definition: TMarker.cxx:117
virtual Color_t GetMarkerColor() const
Return the marker color.
Definition: TAttMarker.h:36
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition: TAttText.h:52
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual TPoint GetBBoxCenter()
Return the center of the BoundingBox as TPoint in pixels.
Definition: TMarker.cxx:355
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0