Jpp 19.3.0-rc.1
the software that should make you happy
Loading...
Searching...
No Matches
JConvexHull2D.cc File Reference

Example program to test convex hull and enclosing circle. More...

#include <string>
#include <iostream>
#include <fstream>
#include <iomanip>
#include <vector>
#include "TROOT.h"
#include "TRandom3.h"
#include "TH2D.h"
#include "TGraph.h"
#include "TMarker.h"
#include "TEllipse.h"
#include "TCanvas.h"
#include "TApplication.h"
#include "JGeometry2D/JVector2D.hh"
#include "JGeometry2D/JCircle2D.hh"
#include "JGeometry2D/JGeometry2DToolkit.hh"
#include "Jeep/JPrint.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Detailed Description

Example program to test convex hull and enclosing circle.

Author
mdejong

Definition in file JConvexHull2D.cc.

Function Documentation

◆ main()

int main ( int argc,
char * argv[] )

Definition at line 32 of file JConvexHull2D.cc.

33{
34 using namespace std;
35
36 string inputFile;
37 string outputFile;
39 ULong_t seed;
40 int debug;
41
42 try {
43
44 JParser<> zap("Example program to test convex hull and enclosing circle.");
45
46 zap['f'] = make_field(inputFile) = "";
47 zap['o'] = make_field(outputFile) = "";
48 zap['n'] = make_field(numberOfPoints) = 0;
49 zap['S'] = make_field(seed) = 0;
50 zap['d'] = make_field(debug) = 0;
51
52 zap(argc, argv);
53 }
54 catch(const exception &error) {
55 FATAL(error.what() << endl);
56 }
57
58 gRandom->SetSeed(seed);
59
60 using namespace JPP;
61
62
63 vector<JVector2D> buffer;
64
65 typedef vector<JVector2D>::const_iterator const_iterator;
67
68
69 if (inputFile != "") {
70
71 ifstream in(inputFile.c_str());
72
73 for (double x, y; in >> x >> y; ) {
74 buffer.push_back(JVector2D(x,y));
75 }
76
77 in.close();
78
79 } else if (numberOfPoints > 0) {
80
81 NOTICE("Seed: " << gRandom->GetSeed() << endl);
82
83 for (int i = 0; i != numberOfPoints; ++i) {
84
85 buffer.push_back(JVector2D(gRandom->Uniform(-1.0, +1.0),
86 gRandom->Uniform(-1.0, +1.0)));
87 }
88
89 if (outputFile != "") {
90
91 ofstream out(outputFile.c_str(), ios::out);
92
93 for (const_iterator i = buffer.begin(); i != buffer.end(); ++i)
94 out << FIXED(7,3) << i->getX() << ' '
95 << FIXED(7,3) << i->getY() << endl;
96
97 out.close();
98 }
99
100 } else {
101
102 FATAL("No points to draw." << endl);
103 }
104
105 for (const_iterator i = buffer.begin(); i != buffer.end(); ++i) {
106 DEBUG(FIXED(7,3) << i->getX() << ' ' << FIXED(7,3) << i->getY() << endl);
107 }
108
109 const JCircle2D circle(buffer.begin(), buffer.end());
110
111 pair<iterator,iterator> hull = getConvexHull2D(buffer.begin(),
112 buffer.end());
113
114
115 DEBUG("upper: "
116 << 0 << ' '
117 << distance(buffer.begin(),hull.first) << endl);
118
119 DEBUG("lower: "
120 << distance(buffer.begin(),hull.first) << ' '
121 << distance(buffer.begin(),hull.second) << endl);
122
123 DEBUG("circle: " << circle.getX() << ' ' << circle.getY() << ' ' << circle.getRadius() << endl);
124
125 double area = getArea2D(buffer.begin(), hull.second);
126
127 NOTICE("Area: " << area << endl);
128
129
130 {
131 int n[] = { 0, 0 };
132
133 for (const_iterator i = buffer.begin(); i != buffer.end(); ++i) {
134 if (!inside2D(buffer.begin(), hull.second, *i)) { ++n[0]; }
135 if (!inside2D(buffer.begin(), hull.first, hull.second, *i)) { ++n[1]; }
136
137 if (!inside2D(buffer.begin(), hull.first, hull.second, *i)) {
138
139 }
140 }
141
142 for (int i = 0; i != sizeof(n)/sizeof(n[0]); ++i) {
143 NOTICE("Number of points outside [" << i << "] " << n[i] << endl);
144 }
145 }
146
147 TApplication* tp = new TApplication("user", NULL, NULL);
148
149 TCanvas cv("cv", "", 400, 400);
150
151 cv.SetFillStyle(4000);
152 cv.SetFillColor(0);
153
154 cv.Divide(1, 1);
155 cv.cd(1);
156
157
158 TEllipse ellipse(circle.getX(), circle.getY(), circle.getRadius());
159
160
161 const Int_t MAX_BUFFER_SIZE = buffer.size() + 1;
162
163 Double_t x[MAX_BUFFER_SIZE];
164 Double_t y[MAX_BUFFER_SIZE];
165
166 int N = 0;
167
168 for (const_iterator i = buffer.begin(); i != buffer.end(); ++i, ++N) {
169 x[N] = i->getX();
170 y[N] = i->getY();
171 }
172
173
174 Double_t xmin = -2.0;
175 Double_t xmax = +2.0;
176 Double_t ymin = -2.0;
177 Double_t ymax = +2.0;
178
179 TH2D h2("h2", "", 1, xmin, xmax, 1, ymin, ymax);
180
181 h2.SetStats(kFALSE);
182 h2.Draw();
183
184
185 ellipse.SetLineWidth(2);
186 ellipse.Draw();
187
188
189 TGraph g(N, x, y);
190
191 g.SetMarkerStyle(20);
192 g.SetMarkerColor(kBlack);
193 g.SetMarkerSize(0.7);
194 g.Draw("P");
195
196
197 // lower hull
198
199 N = distance(buffer.begin(), hull.first);
200
201 TGraph g1(N, x, y);
202
203 g1.SetLineColor(kBlack);
204 g1.SetLineWidth(2);
205 g1.Draw("L");
206
207
208 // upper hull
209
210 const int i = --N;
211
212 N = distance(buffer.begin(), hull.second);
213
214 x[N] = x[0];
215 y[N] = y[0];
216
217 N = N + 1 - i;
218
219 TGraph g2(N, &x[i], &y[i]);
220
221 g2.SetLineColor(kRed);
222 g2.SetLineWidth(2);
223 g2.Draw("L");
224
225
226 cv.Update();
227
228 tp->Run();
229}
string outputFile
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define NOTICE(A)
Definition JMessage.hh:64
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:72
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
Double_t g1(const Double_t x)
Function.
Definition JQuantiles.cc:25
int numberOfPoints
Definition JResultPDF.cc:22
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
Data structure for circle in two dimensions.
Definition JCircle2D.hh:35
Data structure for vector in two dimensions.
Definition JVector2D.hh:34
Utility class to parse command line options.
Definition JParser.hh:1698
const double xmax
const double xmin
double getArea2D(T __begin, T __end)
Get area of a convex polygon.
bool inside2D(T __begin, T __end, const JVector2D &pos)
Check if given point is inside a convex polygon.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
const int n
Definition JPolint.hh:791
Auxiliary data structure for floating point format specification.
Definition JManip.hh:448