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
45{
48 fLastPoint = -1;
49 if (n <= 0) {
50 fN = 0;
51 fLastPoint = -1;
52 fX = fY = nullptr;
53 return;
54 }
55 fN = n;
56 fX = new Double_t [fN];
57 fY = new Double_t [fN];
58}
59
60////////////////////////////////////////////////////////////////////////////////
61/// Constructor.
62
65{
68 fLastPoint = -1;
69 if (n <= 0) {
70 fN = 0;
71 fLastPoint = -1;
72 fX = fY = nullptr;
73 return;
74 }
75 fN = n;
76 fX = new Double_t [fN];
77 fY = new Double_t [fN];
78 if (!x || !y) return;
79 for (Int_t i=0; i<fN;i++) { fX[i] = x[i]; fY[i] = y[i]; }
80 fLastPoint = fN-1;
81}
82
83////////////////////////////////////////////////////////////////////////////////
84/// Constructor.
85
88{
91 fLastPoint = -1;
92 if (n <= 0) {
93 fN = 0;
94 fLastPoint = -1;
95 fX = fY = nullptr;
96 return;
97 }
98 fN = n;
99 fX = new Double_t [fN];
100 fY = new Double_t [fN];
101 if (!x || !y) return;
102 for (Int_t i=0; i<fN;i++) { fX[i] = x[i]; fY[i] = y[i]; }
103 fLastPoint = fN-1;
104}
105
106////////////////////////////////////////////////////////////////////////////////
107///assignment operator
108
110{
111 if(this != &pm)
112 pm.TPolyMarker::Copy(*this);
113
114 return *this;
115}
116
117////////////////////////////////////////////////////////////////////////////////
118/// Destructor.
119
121{
122 if (fX) delete [] fX;
123 if (fY) delete [] fY;
124 fLastPoint = -1;
125}
126
127////////////////////////////////////////////////////////////////////////////////
128/// Copy constructor.
129
130TPolyMarker::TPolyMarker(const TPolyMarker &polymarker) : TObject(polymarker), TAttMarker(polymarker)
131{
132 fN = 0;
133 fX = fY = nullptr;
134 fLastPoint = -1;
135 polymarker.TPolyMarker::Copy(*this);
136}
137
138////////////////////////////////////////////////////////////////////////////////
139// Copy TPolyMarker into provided object
140
142{
143 TObject::Copy(obj);
145 ((TPolyMarker&)obj).fN = fN;
146 // delete first previous existing fX and fY
147 if (((TPolyMarker&)obj).fX) delete [] (((TPolyMarker&)obj).fX);
148 if (((TPolyMarker&)obj).fY) delete [] (((TPolyMarker&)obj).fY);
149 if (fN > 0) {
150 ((TPolyMarker&)obj).fX = new Double_t [fN];
151 ((TPolyMarker&)obj).fY = new Double_t [fN];
152 for (Int_t i=0; i<fN;i++) {
153 ((TPolyMarker&)obj).fX[i] = fX[i];
154 ((TPolyMarker&)obj).fY[i] = fY[i];
155 }
156 } else {
157 ((TPolyMarker&)obj).fX = nullptr;
158 ((TPolyMarker&)obj).fY = nullptr;
159 }
160 ((TPolyMarker&)obj).fOption = fOption;
161 ((TPolyMarker&)obj).fLastPoint = fLastPoint;
162}
163
164////////////////////////////////////////////////////////////////////////////////
165/// Compute distance from point px,py to a polymarker.
166///
167/// Compute the closest distance of approach from point px,py to each point
168/// of the polymarker.
169/// Returns when the distance found is below DistanceMaximum.
170/// The distance is computed in pixels units.
171
173{
174 const Int_t big = 9999;
175
176 // check if point is near one of the points
177 Int_t i, pxp, pyp, d;
178 Int_t distance = big;
179
180 for (i=0;i<Size();i++) {
181 pxp = gPad->XtoAbsPixel(gPad->XtoPad(fX[i]));
182 pyp = gPad->YtoAbsPixel(gPad->YtoPad(fY[i]));
183 d = TMath::Abs(pxp-px) + TMath::Abs(pyp-py);
184 if (d < distance) distance = d;
185 }
186 return distance;
187}
188
189////////////////////////////////////////////////////////////////////////////////
190/// Draw.
191
193{
195}
196
197////////////////////////////////////////////////////////////////////////////////
198/// Draw polymarker.
199
201{
202 TPolyMarker *newpolymarker = new TPolyMarker(n,x,y);
203 TAttMarker::Copy(*newpolymarker);
204 newpolymarker->fOption = fOption;
205 newpolymarker->SetBit(kCanDelete);
206 newpolymarker->AppendPad();
207}
208
209////////////////////////////////////////////////////////////////////////////////
210/// Execute action corresponding to one event.
211///
212/// This member function must be implemented to realize the action
213/// corresponding to the mouse click on the object in the window
214
216{
217}
218
219////////////////////////////////////////////////////////////////////////////////
220/// ls.
221
223{
225 printf("TPolyMarker N=%d\n",fN);
226}
227
228////////////////////////////////////////////////////////////////////////////////
229/// Merge polymarkers in the collection in this polymarker.
230
232{
233 if (!li) return 0;
234 TIter next(li);
235
236 //first loop to count the number of entries
237 TPolyMarker *pm;
238 Int_t npoints = 0;
239 while ((pm = (TPolyMarker*)next())) {
240 if (!pm->InheritsFrom(TPolyMarker::Class())) {
241 Error("Add","Attempt to add object of class: %s to a %s",pm->ClassName(),this->ClassName());
242 return -1;
243 }
244 npoints += pm->Size();
245 }
246
247 //extend this polymarker to hold npoints
248 SetPoint(npoints-1,0,0);
249
250 //merge all polymarkers
251 next.Reset();
252 while ((pm = (TPolyMarker*)next())) {
253 Int_t np = pm->Size();
254 Double_t *x = pm->GetX();
255 Double_t *y = pm->GetY();
256 for (Int_t i=0;i<np;i++) {
257 SetPoint(i,x[i],y[i]);
258 }
259 }
260
261 return npoints;
262}
263
264////////////////////////////////////////////////////////////////////////////////
265/// Paint.
266
268{
270}
271
272////////////////////////////////////////////////////////////////////////////////
273/// Paint polymarker.
274
276{
277 if (n <= 0) return;
278 TAttMarker::Modify(); //Change marker attributes only if necessary
279 Double_t *xx = x;
280 Double_t *yy = y;
281 if (gPad->GetLogx()) {
282 xx = new Double_t[n];
283 for (Int_t ix=0;ix<n;ix++) xx[ix] = gPad->XtoPad(x[ix]);
284 }
285 if (gPad->GetLogy()) {
286 yy = new Double_t[n];
287 for (Int_t iy=0;iy<n;iy++) yy[iy] = gPad->YtoPad(y[iy]);
288 }
289 gPad->PaintPolyMarker(n,xx,yy,option);
290 if (x != xx) delete [] xx;
291 if (y != yy) delete [] yy;
292}
293
294////////////////////////////////////////////////////////////////////////////////
295/// Print polymarker.
296
298{
299 printf("TPolyMarker N=%d\n",fN);
300}
301
302////////////////////////////////////////////////////////////////////////////////
303/// Save primitive as a C++ statement(s) on output stream out.
304
305void TPolyMarker::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
306{
307 char quote = '"';
308 out<<" "<<std::endl;
309 out<<" Double_t *dum = 0;"<<std::endl;
310 if (gROOT->ClassSaved(TPolyMarker::Class())) {
311 out<<" ";
312 } else {
313 out<<" TPolyMarker *";
314 }
315 out<<"pmarker = new TPolyMarker("<<fN<<",dum,dum,"<<quote<<fOption<<quote<<");"<<std::endl;
316
317 SaveMarkerAttributes(out,"pmarker",1,1,1);
318
319 for (Int_t i=0;i<Size();i++) {
320 out<<" pmarker->SetPoint("<<i<<","<<fX[i]<<","<<fY[i]<<");"<<std::endl;
321 }
322 if (!strstr(option, "nodraw")) {
323 out<<" pmarker->Draw("
324 <<quote<<option<<quote<<");"<<std::endl;
325 }
326}
327
328////////////////////////////////////////////////////////////////////////////////
329/// Set point following LastPoint to x, y.
330/// Returns index of the point (new last point).
331
333{
334 fLastPoint++;
336 return fLastPoint;
337}
338
339////////////////////////////////////////////////////////////////////////////////
340/// Set point number n.
341/// if n is greater than the current size, the arrays are automatically
342/// extended
343
345{
346 if (n < 0) return;
347 if (!fX || !fY || n >= fN) {
348 // re-allocate the object
349 Int_t newN = TMath::Max(2*fN,n+1);
350 Double_t *savex = new Double_t [newN];
351 Double_t *savey = new Double_t [newN];
352 if (fX && fN){
353 memcpy(savex,fX,fN*sizeof(Double_t));
354 memset(&savex[fN],0,(newN-fN)*sizeof(Double_t));
355 delete [] fX;
356 }
357 if (fY && fN){
358 memcpy(savey,fY,fN*sizeof(Double_t));
359 memset(&savey[fN],0,(newN-fN)*sizeof(Double_t));
360 delete [] fY;
361 }
362 fX = savex;
363 fY = savey;
364 fN = newN;
365 }
366 fX[n] = x;
367 fY[n] = y;
369}
370
371////////////////////////////////////////////////////////////////////////////////
372/// If n <= 0 the current arrays of points are deleted.
373
375{
376 if (n <= 0) {
377 fN = 0;
378 fLastPoint = -1;
379 delete [] fX;
380 delete [] fY;
381 fX = fY = nullptr;
382 return;
383 }
384 SetPoint(n-1,0,0);
385}
386
387////////////////////////////////////////////////////////////////////////////////
388/// If n <= 0 the current arrays of points are deleted.
389
391{
392 if (n <= 0) {
393 fN = 0;
394 fLastPoint = -1;
395 delete [] fX;
396 delete [] fY;
397 fX = fY = nullptr;
398 return;
399 }
400 fN =n;
401 if (fX) delete [] fX;
402 if (fY) delete [] fY;
403 fX = new Double_t[fN];
404 fY = new Double_t[fN];
405 for (Int_t i=0; i<fN;i++) {
406 if (x) fX[i] = (Double_t)x[i];
407 if (y) fY[i] = (Double_t)y[i];
408 }
409 fOption = option;
410 fLastPoint = fN-1;
411}
412
413////////////////////////////////////////////////////////////////////////////////
414/// If n <= 0 the current arrays of points are deleted.
415
417{
418 if (n <= 0) {
419 fN = 0;
420 fLastPoint = -1;
421 delete [] fX;
422 delete [] fY;
423 fX = fY = nullptr;
424 return;
425 }
426 fN =n;
427 if (fX) delete [] fX;
428 if (fY) delete [] fY;
429 fX = new Double_t[fN];
430 fY = new Double_t[fN];
431 for (Int_t i=0; i<fN;i++) {
432 if (x) fX[i] = x[i];
433 if (y) fY[i] = y[i];
434 }
435 fOption = option;
436 fLastPoint = fN-1;
437}
438
439////////////////////////////////////////////////////////////////////////////////
440/// Stream a class object.
441
443{
444 if (R__b.IsReading()) {
445 UInt_t R__s, R__c;
446 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
447 if (R__v > 1) {
448 R__b.ReadClassBuffer(TPolyMarker::Class(), this, R__v, R__s, R__c);
449 return;
450 }
451 //====process old versions before automatic schema evolution
452 TObject::Streamer(R__b);
454 R__b >> fN;
455 fX = new Double_t[fN];
456 fY = new Double_t[fN];
457 Int_t i;
458 Float_t xold,yold;
459 for (i=0;i<fN;i++) {R__b >> xold; fX[i] = xold;}
460 for (i=0;i<fN;i++) {R__b >> yold; fY[i] = yold;}
461 fOption.Streamer(R__b);
462 R__b.CheckByteCount(R__s, R__c, TPolyMarker::IsA());
463 //====end of old versions
464
465 } else {
467 }
468}
#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:377
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:888
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:207
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:184
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:780
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:525
virtual void Copy(TObject &object) const
Copy this to obj.
Definition TObject.cxx:140
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:987
@ 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:69
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.
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:56
~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:57
virtual void SetPoint(Int_t point, Double_t x, Double_t y)
Set point number n.
TClass * IsA() const override
Definition TPolyMarker.h:71
void ls(Option_t *option="") const override
ls.
static void IndentLevel()
Functions used by ls() to indent an object hierarchy.
Definition TROOT.cxx:2891
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