Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TPolyMarker.cxx
Go to the documentation of this file.
1// @(#)root/hist:$Id$
2// Author: Rene Brun 12/12/94
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 "Riostream.h"
13#include "TROOT.h"
14#include "TBuffer.h"
15#include "TVirtualPad.h"
16#include "TPolyMarker.h"
17#include "TMath.h"
18
20
21
22/** \class TPolyMarker
23 \ingroup Graphs
24A PolyMarker is defined by an array on N points in a 2-D space.
25At each point x[i], y[i] a marker is drawn.
26Marker attributes are managed by TAttMarker.
27See TMarker for the list of possible marker types.
28*/
29
30////////////////////////////////////////////////////////////////////////////////
31/// Default constructor.
32
34{
35 fN = 0;
36 fX = fY = nullptr;
37 fLastPoint = -1;
38}
39
40////////////////////////////////////////////////////////////////////////////////
41/// Constructor.
42
44{
47 fLastPoint = -1;
48 if (n <= 0) {
49 fN = 0;
50 fLastPoint = -1;
51 fX = fY = nullptr;
52 return;
53 }
54 fN = n;
55 fX = new Double_t [fN];
56 fY = new Double_t [fN];
57}
58
59////////////////////////////////////////////////////////////////////////////////
60/// Constructor.
61
63{
66 fLastPoint = -1;
67 if (n <= 0) {
68 fN = 0;
69 fLastPoint = -1;
70 fX = fY = nullptr;
71 return;
72 }
73 fN = n;
74 fX = new Double_t [fN];
75 fY = new Double_t [fN];
76 if (!x || !y) return;
77 for (Int_t i=0; i<fN;i++) { fX[i] = x[i]; fY[i] = y[i]; }
78 fLastPoint = fN-1;
79}
80
81////////////////////////////////////////////////////////////////////////////////
82/// Constructor.
83
85{
88 fLastPoint = -1;
89 if (n <= 0) {
90 fN = 0;
91 fLastPoint = -1;
92 fX = fY = nullptr;
93 return;
94 }
95 fN = n;
96 fX = new Double_t [fN];
97 fY = new Double_t [fN];
98 if (!x || !y) return;
99 for (Int_t i=0; i<fN;i++) { fX[i] = x[i]; fY[i] = y[i]; }
100 fLastPoint = fN-1;
101}
102
103////////////////////////////////////////////////////////////////////////////////
104///assignment operator
105
107{
108 if(this != &pm)
109 pm.TPolyMarker::Copy(*this);
110
111 return *this;
112}
113
114////////////////////////////////////////////////////////////////////////////////
115/// Destructor.
116
118{
119 if (fX) delete [] fX;
120 if (fY) delete [] fY;
121 fLastPoint = -1;
122}
123
124////////////////////////////////////////////////////////////////////////////////
125/// Copy constructor.
126
127TPolyMarker::TPolyMarker(const TPolyMarker &polymarker) : TObject(polymarker), TAttMarker(polymarker)
128{
129 fN = 0;
130 fX = fY = nullptr;
131 fLastPoint = -1;
132 polymarker.TPolyMarker::Copy(*this);
133}
134
135////////////////////////////////////////////////////////////////////////////////
136// Copy TPolyMarker into provided object
137
139{
140 TObject::Copy(obj);
142 ((TPolyMarker&)obj).fN = fN;
143 // delete first previous existing fX and fY
144 if (((TPolyMarker&)obj).fX) delete [] (((TPolyMarker&)obj).fX);
145 if (((TPolyMarker&)obj).fY) delete [] (((TPolyMarker&)obj).fY);
146 if (fN > 0) {
147 ((TPolyMarker&)obj).fX = new Double_t [fN];
148 ((TPolyMarker&)obj).fY = new Double_t [fN];
149 for (Int_t i=0; i<fN;i++) {
150 ((TPolyMarker&)obj).fX[i] = fX[i];
151 ((TPolyMarker&)obj).fY[i] = fY[i];
152 }
153 } else {
154 ((TPolyMarker&)obj).fX = nullptr;
155 ((TPolyMarker&)obj).fY = nullptr;
156 }
157 ((TPolyMarker&)obj).fOption = fOption;
158 ((TPolyMarker&)obj).fLastPoint = fLastPoint;
159}
160
161////////////////////////////////////////////////////////////////////////////////
162/// Compute distance from point px,py to a polymarker.
163///
164/// Compute the closest distance of approach from point px,py to each point
165/// of the polymarker.
166/// Returns when the distance found is below DistanceMaximum.
167/// The distance is computed in pixels units.
168
170{
171 const Int_t big = 9999;
172
173 // check if point is near one of the points
174 Int_t i, pxp, pyp, d;
175 Int_t distance = big;
176
177 for (i=0;i<Size();i++) {
178 pxp = gPad->XtoAbsPixel(gPad->XtoPad(fX[i]));
179 pyp = gPad->YtoAbsPixel(gPad->YtoPad(fY[i]));
180 d = TMath::Abs(pxp-px) + TMath::Abs(pyp-py);
181 if (d < distance) distance = d;
182 }
183 return distance;
184}
185
186////////////////////////////////////////////////////////////////////////////////
187/// Draw.
188
190{
192}
193
194////////////////////////////////////////////////////////////////////////////////
195/// Draw polymarker.
196
198{
199 TPolyMarker *newpolymarker = new TPolyMarker(n,x,y);
200 TAttMarker::Copy(*newpolymarker);
201 newpolymarker->fOption = fOption;
202 newpolymarker->SetBit(kCanDelete);
203 newpolymarker->AppendPad();
204}
205
206////////////////////////////////////////////////////////////////////////////////
207/// Execute action corresponding to one event.
208///
209/// This member function must be implemented to realize the action
210/// corresponding to the mouse click on the object in the window
211
213{
214}
215
216////////////////////////////////////////////////////////////////////////////////
217/// ls.
218
220{
222 printf("TPolyMarker N=%d\n",fN);
223}
224
225////////////////////////////////////////////////////////////////////////////////
226/// Merge polymarkers in the collection in this polymarker.
227
229{
230 if (!li) return 0;
231 TIter next(li);
232
233 //first loop to count the number of entries
234 TPolyMarker *pm;
235 Int_t npoints = 0;
236 while ((pm = (TPolyMarker*)next())) {
237 if (!pm->InheritsFrom(TPolyMarker::Class())) {
238 Error("Add","Attempt to add object of class: %s to a %s",pm->ClassName(),this->ClassName());
239 return -1;
240 }
241 npoints += pm->Size();
242 }
243
244 //extend this polymarker to hold npoints
245 SetPoint(npoints-1,0,0);
246
247 //merge all polymarkers
248 next.Reset();
249 while ((pm = (TPolyMarker*)next())) {
250 Int_t np = pm->Size();
251 Double_t *x = pm->GetX();
252 Double_t *y = pm->GetY();
253 for (Int_t i=0;i<np;i++) {
254 SetPoint(i,x[i],y[i]);
255 }
256 }
257
258 return npoints;
259}
260
261////////////////////////////////////////////////////////////////////////////////
262/// Paint.
263
265{
267}
268
269////////////////////////////////////////////////////////////////////////////////
270/// Paint polymarker.
271
273{
274 if (n <= 0) return;
275 TAttMarker::Modify(); //Change marker attributes only if necessary
276 Double_t *xx = x;
277 Double_t *yy = y;
278 if (gPad->GetLogx()) {
279 xx = new Double_t[n];
280 for (Int_t ix=0;ix<n;ix++) xx[ix] = gPad->XtoPad(x[ix]);
281 }
282 if (gPad->GetLogy()) {
283 yy = new Double_t[n];
284 for (Int_t iy=0;iy<n;iy++) yy[iy] = gPad->YtoPad(y[iy]);
285 }
286 gPad->PaintPolyMarker(n,xx,yy,option);
287 if (x != xx) delete [] xx;
288 if (y != yy) delete [] yy;
289}
290
291////////////////////////////////////////////////////////////////////////////////
292/// Print polymarker.
293
295{
296 printf("TPolyMarker N=%d\n",fN);
297}
298
299////////////////////////////////////////////////////////////////////////////////
300/// Save primitive as a C++ statement(s) on output stream out.
301
302void TPolyMarker::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
303{
304 char quote = '"';
305 out<<" "<<std::endl;
306 out<<" Double_t *dum = 0;"<<std::endl;
307 if (gROOT->ClassSaved(TPolyMarker::Class())) {
308 out<<" ";
309 } else {
310 out<<" TPolyMarker *";
311 }
312 out<<"pmarker = new TPolyMarker("<<fN<<",dum,dum,"<<quote<<fOption<<quote<<");"<<std::endl;
313
314 SaveMarkerAttributes(out,"pmarker",1,1,1);
315
316 for (Int_t i=0;i<Size();i++) {
317 out<<" pmarker->SetPoint("<<i<<","<<fX[i]<<","<<fY[i]<<");"<<std::endl;
318 }
319 if (!strstr(option, "nodraw")) {
320 out<<" pmarker->Draw("
321 <<quote<<option<<quote<<");"<<std::endl;
322 }
323}
324
325////////////////////////////////////////////////////////////////////////////////
326/// Set point following LastPoint to x, y.
327/// Returns index of the point (new last point).
328
330{
331 fLastPoint++;
333 return fLastPoint;
334}
335
336////////////////////////////////////////////////////////////////////////////////
337/// Set point number n.
338/// if n is greater than the current size, the arrays are automatically
339/// extended
340
342{
343 if (n < 0) return;
344 if (!fX || !fY || n >= fN) {
345 // re-allocate the object
346 Int_t newN = TMath::Max(2*fN,n+1);
347 Double_t *savex = new Double_t [newN];
348 Double_t *savey = new Double_t [newN];
349 if (fX && fN){
350 memcpy(savex,fX,fN*sizeof(Double_t));
351 memset(&savex[fN],0,(newN-fN)*sizeof(Double_t));
352 delete [] fX;
353 }
354 if (fY && fN){
355 memcpy(savey,fY,fN*sizeof(Double_t));
356 memset(&savey[fN],0,(newN-fN)*sizeof(Double_t));
357 delete [] fY;
358 }
359 fX = savex;
360 fY = savey;
361 fN = newN;
362 }
363 fX[n] = x;
364 fY[n] = y;
366}
367
368////////////////////////////////////////////////////////////////////////////////
369/// If n <= 0 the current arrays of points are deleted.
370
372{
373 if (n <= 0) {
374 fN = 0;
375 fLastPoint = -1;
376 delete [] fX;
377 delete [] fY;
378 fX = fY = nullptr;
379 return;
380 }
381 SetPoint(n-1,0,0);
382}
383
384////////////////////////////////////////////////////////////////////////////////
385/// If n <= 0 the current arrays of points are deleted.
386
388{
389 if (n <= 0) {
390 fN = 0;
391 fLastPoint = -1;
392 delete [] fX;
393 delete [] fY;
394 fX = fY = nullptr;
395 return;
396 }
397 fN =n;
398 if (fX) delete [] fX;
399 if (fY) delete [] fY;
400 fX = new Double_t[fN];
401 fY = new Double_t[fN];
402 for (Int_t i=0; i<fN;i++) {
403 if (x) fX[i] = (Double_t)x[i];
404 if (y) fY[i] = (Double_t)y[i];
405 }
406 fOption = option;
407 fLastPoint = fN-1;
408}
409
410////////////////////////////////////////////////////////////////////////////////
411/// If n <= 0 the current arrays of points are deleted.
412
414{
415 if (n <= 0) {
416 fN = 0;
417 fLastPoint = -1;
418 delete [] fX;
419 delete [] fY;
420 fX = fY = nullptr;
421 return;
422 }
423 fN =n;
424 if (fX) delete [] fX;
425 if (fY) delete [] fY;
426 fX = new Double_t[fN];
427 fY = new Double_t[fN];
428 for (Int_t i=0; i<fN;i++) {
429 if (x) fX[i] = x[i];
430 if (y) fY[i] = y[i];
431 }
432 fOption = option;
433 fLastPoint = fN-1;
434}
435
436////////////////////////////////////////////////////////////////////////////////
437/// Stream a class object.
438
440{
441 if (R__b.IsReading()) {
442 UInt_t R__s, R__c;
443 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
444 if (R__v > 1) {
445 R__b.ReadClassBuffer(TPolyMarker::Class(), this, R__v, R__s, R__c);
446 return;
447 }
448 //====process old versions before automatic schema evolution
449 TObject::Streamer(R__b);
451 R__b >> fN;
452 fX = new Double_t[fN];
453 fY = new Double_t[fN];
454 Int_t i;
455 Float_t xold,yold;
456 for (i=0;i<fN;i++) {R__b >> xold; fX[i] = xold;}
457 for (i=0;i<fN;i++) {R__b >> yold; fY[i] = yold;}
458 fOption.Streamer(R__b);
459 R__b.CheckByteCount(R__s, R__c, TPolyMarker::IsA());
460 //====end of old versions
461
462 } else {
464 }
465}
#define d(i)
Definition RSha256.hxx:102
short Version_t
Definition RtypesCore.h:65
float Float_t
Definition RtypesCore.h:57
double Double_t
Definition RtypesCore.h:59
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:382
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np
#define gROOT
Definition TROOT.h:406
#define gPad
Marker Attributes class.
Definition TAttMarker.h:19
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.
virtual void Modify()
Change current marker attributes if necessary.
void Copy(TAttMarker &attmarker) const
Copy this marker attributes to a new TAttMarker.
virtual void Streamer(TBuffer &)
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual Version_t ReadVersion(UInt_t *start=nullptr, UInt_t *bcnt=nullptr, const TClass *cl=nullptr)=0
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=nullptr)=0
Bool_t IsReading() const
Definition TBuffer.h:86
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
Collection abstract base class.
Definition TCollection.h:65
void Reset()
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Streamer(TBuffer &)
Stream an object of class TObject.
Definition TObject.cxx:906
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:225
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:202
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:798
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:542
virtual void Copy(TObject &object) const
Copy this to obj.
Definition TObject.cxx:158
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1005
@ kCanDelete
if object in a list can be deleted
Definition TObject.h:62
A PolyMarker is defined by an array on N points in a 2-D space.
Definition TPolyMarker.h:31
TPolyMarker & operator=(const TPolyMarker &)
assignment operator
virtual void PaintPolyMarker(Int_t n, Double_t *x, Double_t *y, Option_t *option="")
Paint polymarker.
void Copy(TObject &polymarker) const override
Copy this to obj.
void Print(Option_t *option="") const override
Print polymarker.
Double_t * fY
[fN] Array of Y coordinates
Definition TPolyMarker.h:36
virtual Int_t Size() const
Definition TPolyMarker.h:76
virtual Int_t Merge(TCollection *list)
Merge polymarkers in the collection in this polymarker.
TPolyMarker()
Default constructor.
virtual void DrawPolyMarker(Int_t n, Double_t *x, Double_t *y, Option_t *option="")
Draw polymarker.
virtual void SetPolyMarker(Int_t n)
If n <= 0 the current arrays of points are deleted.
void Streamer(TBuffer &) override
Stream a class object.
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitive as a C++ statement(s) on output stream out.
Int_t fN
Number of points internally reserved (not necessarily used)
Definition TPolyMarker.h:33
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
Double_t * GetX() const
Definition TPolyMarker.h:63
~TPolyMarker() override
Destructor.
virtual Int_t SetNextPoint(Double_t x, Double_t y)
Set point following LastPoint to x, y.
static TClass * Class()
Int_t fLastPoint
The index of the last filled point.
Definition TPolyMarker.h:34
TString fOption
Options.
Definition TPolyMarker.h:37
void Paint(Option_t *option="") override
Paint.
Double_t * fX
[fN] Array of X coordinates
Definition TPolyMarker.h:35
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Compute distance from point px,py to a polymarker.
void Draw(Option_t *option="") override
Draw.
Double_t * GetY() const
Definition TPolyMarker.h:64
virtual void SetPoint(Int_t point, Double_t x, Double_t y)
Set point number n.
TClass * IsA() const override
Definition TPolyMarker.h:78
void ls(Option_t *option="") const override
ls.
static void IndentLevel()
Functions used by ls() to indent an object hierarchy.
Definition TROOT.cxx:2895
virtual void Streamer(TBuffer &)
Stream a string object.
Definition TString.cxx:1412
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123