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
19
20
21/** \class TPolyMarker
22 \ingroup Graphs
23A PolyMarker is defined by an array on N points in a 2-D space.
24At each point x[i], y[i] a marker is drawn.
25Marker attributes are managed by TAttMarker.
26See TMarker for the list of possible marker types.
27*/
28
29////////////////////////////////////////////////////////////////////////////////
30/// Default constructor.
31
33{
34 fN = 0;
35 fX = fY = nullptr;
36 fLastPoint = -1;
37}
38
39////////////////////////////////////////////////////////////////////////////////
40/// Constructor.
41
43{
46 fLastPoint = -1;
47 if (n <= 0) {
48 fN = 0;
49 fLastPoint = -1;
50 fX = fY = nullptr;
51 return;
52 }
53 fN = n;
54 fX = new Double_t [fN];
55 fY = new Double_t [fN];
56}
57
58////////////////////////////////////////////////////////////////////////////////
59/// Constructor.
60
62{
65 fLastPoint = -1;
66 if (n <= 0) {
67 fN = 0;
68 fLastPoint = -1;
69 fX = fY = nullptr;
70 return;
71 }
72 fN = n;
73 fX = new Double_t [fN];
74 fY = new Double_t [fN];
75 if (!x || !y) return;
76 for (Int_t i=0; i<fN;i++) { fX[i] = x[i]; fY[i] = y[i]; }
77 fLastPoint = fN-1;
78}
79
80////////////////////////////////////////////////////////////////////////////////
81/// Constructor.
82
84{
87 fLastPoint = -1;
88 if (n <= 0) {
89 fN = 0;
90 fLastPoint = -1;
91 fX = fY = nullptr;
92 return;
93 }
94 fN = n;
95 fX = new Double_t [fN];
96 fY = new Double_t [fN];
97 if (!x || !y) return;
98 for (Int_t i=0; i<fN;i++) { fX[i] = x[i]; fY[i] = y[i]; }
99 fLastPoint = fN-1;
100}
101
102////////////////////////////////////////////////////////////////////////////////
103///assignment operator
104
106{
107 if(this != &pm)
108 pm.TPolyMarker::Copy(*this);
109
110 return *this;
111}
112
113////////////////////////////////////////////////////////////////////////////////
114/// Destructor.
115
117{
118 if (fX) delete [] fX;
119 if (fY) delete [] fY;
120 fLastPoint = -1;
121}
122
123////////////////////////////////////////////////////////////////////////////////
124/// Copy constructor.
125
127{
128 fN = 0;
129 fX = fY = nullptr;
130 fLastPoint = -1;
131 polymarker.TPolyMarker::Copy(*this);
132}
133
134////////////////////////////////////////////////////////////////////////////////
135// Copy TPolyMarker into provided object
136
138{
139 TObject::Copy(obj);
141 ((TPolyMarker&)obj).fN = fN;
142 // delete first previous existing fX and fY
143 if (((TPolyMarker&)obj).fX) delete [] (((TPolyMarker&)obj).fX);
144 if (((TPolyMarker&)obj).fY) delete [] (((TPolyMarker&)obj).fY);
145 if (fN > 0) {
146 ((TPolyMarker&)obj).fX = new Double_t [fN];
147 ((TPolyMarker&)obj).fY = new Double_t [fN];
148 for (Int_t i=0; i<fN;i++) {
149 ((TPolyMarker&)obj).fX[i] = fX[i];
150 ((TPolyMarker&)obj).fY[i] = fY[i];
151 }
152 } else {
153 ((TPolyMarker&)obj).fX = nullptr;
154 ((TPolyMarker&)obj).fY = nullptr;
155 }
156 ((TPolyMarker&)obj).fOption = fOption;
157 ((TPolyMarker&)obj).fLastPoint = fLastPoint;
158}
159
160////////////////////////////////////////////////////////////////////////////////
161/// Compute distance from point px,py to a polymarker.
162///
163/// Compute the closest distance of approach from point px,py to each point
164/// of the polymarker.
165/// Returns when the distance found is below DistanceMaximum.
166/// The distance is computed in pixels units.
167
169{
170 const Int_t big = 9999;
171
172 // check if point is near one of the points
173 Int_t i, pxp, pyp, d;
175
176 for (i=0;i<Size();i++) {
177 pxp = gPad->XtoAbsPixel(gPad->XtoPad(fX[i]));
178 pyp = gPad->YtoAbsPixel(gPad->YtoPad(fY[i]));
179 d = TMath::Abs(pxp-px) + TMath::Abs(pyp-py);
180 if (d < distance) distance = d;
181 }
182 return distance;
183}
184
185////////////////////////////////////////////////////////////////////////////////
186/// Draw.
187
192
193////////////////////////////////////////////////////////////////////////////////
194/// Draw polymarker.
195
204
205////////////////////////////////////////////////////////////////////////////////
206/// Execute action corresponding to one event.
207///
208/// This member function must be implemented to realize the action
209/// corresponding to the mouse click on the object in the window
210
214
215////////////////////////////////////////////////////////////////////////////////
216/// ls.
217
219{
221 printf("TPolyMarker N=%d\n",fN);
222}
223
224////////////////////////////////////////////////////////////////////////////////
225/// Merge polymarkers in the collection in this polymarker.
226
228{
229 if (!li) return 0;
230 TIter next(li);
231
232 //first loop to count the number of entries
234 Int_t npoints = 0;
235 while ((pm = (TPolyMarker*)next())) {
236 if (!pm->InheritsFrom(TPolyMarker::Class())) {
237 Error("Add","Attempt to add object of class: %s to a %s",pm->ClassName(),this->ClassName());
238 return -1;
239 }
240 npoints += pm->Size();
241 }
242
243 //extend this polymarker to hold npoints
244 SetPoint(npoints-1,0,0);
245
246 //merge all polymarkers
247 next.Reset();
248 while ((pm = (TPolyMarker*)next())) {
249 Int_t np = pm->Size();
250 Double_t *x = pm->GetX();
251 Double_t *y = pm->GetY();
252 for (Int_t i=0;i<np;i++) {
253 SetPoint(i,x[i],y[i]);
254 }
255 }
256
257 return npoints;
258}
259
260////////////////////////////////////////////////////////////////////////////////
261/// Paint.
262
267
268////////////////////////////////////////////////////////////////////////////////
269/// Paint polymarker.
270
272{
273 if (n <= 0) return;
274 TAttMarker::Modify(); //Change marker attributes only if necessary
275 Double_t *xx = x;
276 Double_t *yy = y;
277 if (gPad->GetLogx()) {
278 xx = new Double_t[n];
279 for (Int_t ix=0;ix<n;ix++) xx[ix] = gPad->XtoPad(x[ix]);
280 }
281 if (gPad->GetLogy()) {
282 yy = new Double_t[n];
283 for (Int_t iy=0;iy<n;iy++) yy[iy] = gPad->YtoPad(y[iy]);
284 }
285 gPad->PaintPolyMarker(n,xx,yy,option);
286 if (x != xx) delete [] xx;
287 if (y != yy) delete [] yy;
288}
289
290////////////////////////////////////////////////////////////////////////////////
291/// Print polymarker.
292
294{
295 printf("TPolyMarker N=%d\n",fN);
296}
297
298////////////////////////////////////////////////////////////////////////////////
299/// Save primitive as a C++ statement(s) on output stream out.
300
302{
303 TString args;
304 if (Size() > 0) {
305 TString arrxname = SavePrimitiveVector(out, "pmarker", Size(), fX, kTRUE);
306 TString arryname = SavePrimitiveVector(out, "pmarker", Size(), fY);
307 args.Form("%d, %s.data(), %s.data(), \"", Size(), arrxname.Data(), arryname.Data());
308 } else
309 args = "0, \"";
310 args.Append(TString(fOption).ReplaceSpecialCppChars() + "\"");
311
312 SavePrimitiveConstructor(out, Class(), "pmarker", args, Size() == 0);
313
314 SaveMarkerAttributes(out, "pmarker", 1, 1, 1);
315
316 SavePrimitiveDraw(out, "pmarker", option);
317}
318
319////////////////////////////////////////////////////////////////////////////////
320/// Set point following LastPoint to x, y.
321/// Returns index of the point (new last point).
322
329
330////////////////////////////////////////////////////////////////////////////////
331/// Set point number n.
332/// if n is greater than the current size, the arrays are automatically
333/// extended
334
336{
337 if (n < 0) return;
338 if (!fX || !fY || n >= fN) {
339 // re-allocate the object
340 Int_t newN = TMath::Max(2*fN,n+1);
341 Double_t *savex = new Double_t [newN];
342 Double_t *savey = new Double_t [newN];
343 if (fX && fN){
344 memcpy(savex,fX,fN*sizeof(Double_t));
345 memset(&savex[fN],0,(newN-fN)*sizeof(Double_t));
346 delete [] fX;
347 }
348 if (fY && fN){
349 memcpy(savey,fY,fN*sizeof(Double_t));
350 memset(&savey[fN],0,(newN-fN)*sizeof(Double_t));
351 delete [] fY;
352 }
353 fX = savex;
354 fY = savey;
355 fN = newN;
356 }
357 fX[n] = x;
358 fY[n] = y;
360}
361
362////////////////////////////////////////////////////////////////////////////////
363/// If n <= 0 the current arrays of points are deleted.
364
366{
367 if (n <= 0) {
368 fN = 0;
369 fLastPoint = -1;
370 delete [] fX;
371 delete [] fY;
372 fX = fY = nullptr;
373 return;
374 }
375 SetPoint(n-1,0,0);
376}
377
378////////////////////////////////////////////////////////////////////////////////
379/// If n <= 0 the current arrays of points are deleted.
380
382{
383 if (n <= 0) {
384 fN = 0;
385 fLastPoint = -1;
386 delete [] fX;
387 delete [] fY;
388 fX = fY = nullptr;
389 return;
390 }
391 fN =n;
392 if (fX) delete [] fX;
393 if (fY) delete [] fY;
394 fX = new Double_t[fN];
395 fY = new Double_t[fN];
396 for (Int_t i=0; i<fN;i++) {
397 if (x) fX[i] = (Double_t)x[i];
398 if (y) fY[i] = (Double_t)y[i];
399 }
400 fOption = option;
401 fLastPoint = fN-1;
402}
403
404////////////////////////////////////////////////////////////////////////////////
405/// If n <= 0 the current arrays of points are deleted.
406
408{
409 if (n <= 0) {
410 fN = 0;
411 fLastPoint = -1;
412 delete [] fX;
413 delete [] fY;
414 fX = fY = nullptr;
415 return;
416 }
417 fN =n;
418 if (fX) delete [] fX;
419 if (fY) delete [] fY;
420 fX = new Double_t[fN];
421 fY = new Double_t[fN];
422 for (Int_t i=0; i<fN;i++) {
423 if (x) fX[i] = x[i];
424 if (y) fY[i] = y[i];
425 }
426 fOption = option;
427 fLastPoint = fN-1;
428}
429
430////////////////////////////////////////////////////////////////////////////////
431/// Stream a class object.
432
434{
435 if (R__b.IsReading()) {
437 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
438 if (R__v > 1) {
439 R__b.ReadClassBuffer(TPolyMarker::Class(), this, R__v, R__s, R__c);
440 return;
441 }
442 //====process old versions before automatic schema evolution
445 R__b >> fN;
446 fX = new Double_t[fN];
447 fY = new Double_t[fN];
448 Int_t i;
450 for (i=0;i<fN;i++) {R__b >> xold; fX[i] = xold;}
451 for (i=0;i<fN;i++) {R__b >> yold; fY[i] = yold;}
453 R__b.CheckByteCount(R__s, R__c, TPolyMarker::IsA());
454 //====end of old versions
455
456 } else {
457 R__b.WriteClassBuffer(TPolyMarker::Class(),this);
458 }
459}
#define d(i)
Definition RSha256.hxx:102
short Version_t
Class version identifier (short)
Definition RtypesCore.h:79
float Float_t
Float 4 bytes (float)
Definition RtypesCore.h:71
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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 gPad
Marker Attributes class.
Definition TAttMarker.h:20
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
Collection abstract base class.
Definition TCollection.h:65
void Reset()
Mother of all ROOT objects.
Definition TObject.h:41
static TString SavePrimitiveVector(std::ostream &out, const char *prefix, Int_t len, Double_t *arr, Bool_t empty_line=kFALSE)
Save array in the output stream "out" as vector.
Definition TObject.cxx:788
virtual void Streamer(TBuffer &)
Stream an object of class TObject.
Definition TObject.cxx:972
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:203
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:864
virtual void Copy(TObject &object) const
Copy this to obj.
Definition TObject.cxx:159
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1071
static void SavePrimitiveDraw(std::ostream &out, const char *variable_name, Option_t *option=nullptr)
Save invocation of primitive Draw() method Skipped if option contains "nodraw" string.
Definition TObject.cxx:822
static void SavePrimitiveConstructor(std::ostream &out, TClass *cl, const char *variable_name, const char *constructor_agrs="", Bool_t empty_line=kTRUE)
Save object constructor in the output stream "out".
Definition TObject.cxx:771
@ kCanDelete
if object in a list can be deleted
Definition TObject.h:68
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.
~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.
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:2898
Basic string class.
Definition TString.h:138
virtual void Streamer(TBuffer &)
Stream a string object.
Definition TString.cxx:1418
TString & Append(const char *cs)
Definition TString.h:580
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition TString.cxx:2362
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:251
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:124