25 void th2polyEurope(Int_t npoints=500000)
32 Double_t R = (lat2-lat1)/(lon2-lon1);
34 Int_t H = (Int_t)(R*800);
35 gStyle->SetStatX(0.28);
36 gStyle->SetStatY(0.45);
37 gStyle->SetStatW(0.15);
40 TCanvas *ce =
new TCanvas(
"ce",
"ce",0,0,W,H);
41 ce->ToggleEventStatus();
48 const char *countries[nx] = {
49 "france",
"spain",
"sweden",
"germany",
"finland",
50 "norway",
"poland",
"italy",
"yugoslavia",
"united_kingdom",
51 "romania",
"belarus",
"greece",
"czechoslovakia",
"bulgaria",
52 "iceland",
"hungary",
"portugal",
"austria",
"ireland",
53 "lithuania",
"latvia",
"estonia",
"denmark",
"netherlands",
54 "switzerland",
"moldova",
"belgium",
"albania",
"cyprus",
55 "luxembourg",
"andorra",
"malta",
"liechtenstein",
"san_marino",
"monaco" };
56 Float_t surfaces[nx] = {
57 547030, 505580, 449964, 357021, 338145,
58 324220, 312685, 301230, 255438, 244820,
59 237500, 207600, 131940, 127711, 110910,
60 103000, 93030, 89242, 83870, 70280,
61 65200, 64589, 45226, 43094, 41526,
62 41290, 33843, 30528, 28748, 9250,
63 2586, 468, 316, 160, 61, 2};
65 TH1F *h =
new TH1F(
"h",
"Countries surfaces (in km^{2})",3,0,3);
66 for (i=0; i<nx; i++) h->Fill(countries[i], surfaces[i]);
69 TFile::SetCacheFileDir(
".");
71 f = TFile::Open(
"http://root.cern.ch/files/europe.root",
"cacheread");
74 printf(
"Cannot access europe.root. Is internet working ?\n");
78 TH2Poly *p =
new TH2Poly(
80 "Europe (bin contents are normalized to the surfaces in km^{2})",
82 p->GetXaxis()->SetNdivisions(520);
83 p->GetXaxis()->SetTitle(
"longitude");
84 p->GetYaxis()->SetTitle(
"latitude");
90 TIter nextkey(gDirectory->GetListOfKeys());
91 while ((key = (TKey*)nextkey())) {
92 TObject *obj = key->ReadObj();
93 if (obj->InheritsFrom(
"TMultiGraph")) {
94 mg = (TMultiGraph*)obj;
100 Double_t longitude, latitude;
101 Double_t x, y, pi4 = TMath::Pi()/4, alpha = TMath::Pi()/360;
103 gBenchmark->Start(
"Partitioning");
104 p->ChangePartition(100, 100);
105 gBenchmark->Show(
"Partitioning");
108 gBenchmark->Start(
"Filling");
109 for (i=0; i<npoints; i++) {
110 longitude = r.Uniform(lon1,lon2);
111 latitude = r.Uniform(lat1,lat2);
113 y = 38*TMath::Log(TMath::Tan(pi4+alpha*latitude));
116 gBenchmark->Show(
"Filling");
118 Int_t nbins = p->GetNumberOfBins();
119 Double_t maximum = p->GetMaximum();
123 TH1F *h2 = (TH1F *)h->Clone(
"h2");
125 for (j=0; j<nx; j++) {
126 for (i=0; i<nbins; i++) {
127 if (strstr(countries[j],p->GetBinName(i+1))) {
128 h2->Fill(countries[j],p->GetBinContent(i+1));
129 h2->SetBinError(j, p->GetBinError(i+1));
135 Double_t scale = surfaces[0]/maximum;
136 for (i=0; i<nbins; i++) p->SetBinContent(i+1, scale*p->GetBinContent(i+1));
138 gStyle->SetOptStat(1111);
141 TCanvas *c1 =
new TCanvas(
"c1",
"c1",W+10,0,W-20,H);
142 c1->SetRightMargin(0.047);
144 scale = h->GetMaximum()/h2->GetMaximum();
147 h->SetLineColor(kRed-3);
149 h->SetMarkerStyle(20);
150 h->SetMarkerColor(kBlue);
151 h->SetMarkerSize(0.8);
153 h->GetXaxis()->SetLabelFont(42);
154 h->GetXaxis()->SetLabelSize(0.03);
155 h->GetYaxis()->SetLabelFont(42);
158 Double_t scale2=TMath::Sqrt(scale);
159 for (i=0; i<nx; i++) h2->SetBinError(i+1, scale2*h2->GetBinError(i+1));
161 h2->SetMarkerStyle(20);
162 h2->SetMarkerSize(0.8);
164 TLegend *leg =
new TLegend(0.5,0.67,0.92,0.8,NULL,
"NDC");
165 leg->SetTextFont(42);
166 leg->SetTextSize(0.025);
167 leg->AddEntry(h,
"Real countries surfaces from Wikipedia (in km^{2})",
"lp");
168 leg->AddEntry(h2,
"Countries surfaces from TH2Poly (with errors)",
"lp");
172 Double_t wikiSum = h->Integral();
173 Double_t polySum = h2->Integral();
174 Double_t error = TMath::Abs(wikiSum-polySum)/wikiSum;
175 printf(
"THPoly Europe surface estimation error wrt wikipedia = %f per cent when using %d points\n",100*error,npoints);