Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoChecker.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Andrei Gheata 01/11/01
3// CheckGeometry(), CheckOverlaps() by Mihaela Gheata
4
5/*************************************************************************
6 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
7 * All rights reserved. *
8 * *
9 * For the licensing terms see $ROOTSYS/LICENSE. *
10 * For the list of contributors see $ROOTSYS/README/CREDITS. *
11 *************************************************************************/
12
13/** \class TGeoChecker
14\ingroup Geometry_painter
15
16Geometry checking package.
17
18TGeoChecker class provides several geometry checking methods. There are two
19types of tests that can be performed. One is based on random sampling or
20ray-tracing and provides a visual check on how navigation methods work for
21a given geometry. The second actually checks the validity of the geometry
22definition in terms of overlapping/extruding objects. Both types of checks
23can be done for a given branch (starting with a given volume) as well as for
24the geometry as a whole.
25
26#### TGeoChecker::CheckPoint(Double_t x, Double_t y, Double_t z)
27
28This method can be called directly from the TGeoManager class and print a
29report on how a given point is classified by the modeller: which is the
30full path to the deepest node containing it, and the (under)estimation
31of the distance to the closest boundary (safety).
32
33#### TGeoChecker::RandomPoints(Int_t npoints)
34
35Can be called from TGeoVolume class. It first draws the volume and its
36content with the current visualization settings. Then randomly samples points
37in its bounding box, plotting in the geometry display only the points
38classified as belonging to visible volumes.
39
40#### TGeoChecker::RandomRays(Int_t nrays, Double_t startx, starty, startz)
41
42Can be called and acts in the same way as the previous, but instead of points,
43rays having random isotropic directions are generated from the given point.
44A raytracing algorithm propagates all rays until they exit geometry, plotting
45all segments crossing visible nodes in the same color as these.
46
47#### TGeoChecker::Test(Int_t npoints)
48
49Implementation of TGeoManager::Test(). Computes the time for the modeller
50to find out "Where am I?" for a given number of random points.
51
52#### TGeoChecker::LegoPlot(ntheta, themin, themax, nphi, phimin, phimax,...)
53
54Implementation of TGeoVolume::LegoPlot(). Draws a spherical radiation length
55lego plot for a given volume, in a given theta/phi range.
56
57#### TGeoChecker::Weight(Double_t precision)
58
59Implementation of TGeoVolume::Weight(). Estimates the total weight of a given
60volume by material sampling. Accepts as input the desired precision.
61*/
62
63#include "TVirtualPad.h"
64#include "TCanvas.h"
65#include "TStyle.h"
66#include "TFile.h"
67#include "TNtuple.h"
68#include "TH2.h"
69#include "TRandom3.h"
70#include "TPolyMarker3D.h"
71#include "TPolyLine3D.h"
72#include "TStopwatch.h"
73
74#include "TGeoVoxelFinder.h"
75#include "TGeoBBox.h"
76#include "TGeoPcon.h"
77#include "TGeoTessellated.h"
78#include "TGeoManager.h"
79#include "TGeoOverlap.h"
80#include "TGeoPainter.h"
81#include "TGeoChecker.h"
82
83#include "TBuffer3D.h"
84#include "TBuffer3DTypes.h"
85#include "TMath.h"
86
87#include <cstdlib>
88
89// statics and globals
90
92
93////////////////////////////////////////////////////////////////////////////////
94/// Default constructor
95
97 :TObject(),
98 fGeoManager(NULL),
99 fVsafe(NULL),
100 fBuff1(NULL),
101 fBuff2(NULL),
102 fFullCheck(kFALSE),
103 fVal1(NULL),
104 fVal2(NULL),
105 fFlags(NULL),
106 fTimer(NULL),
107 fSelectedNode(NULL),
108 fNchecks(0),
109 fNmeshPoints(1000)
110{
111}
112
113////////////////////////////////////////////////////////////////////////////////
114/// Constructor for a given geometry
115
117 :TObject(),
118 fGeoManager(geom),
119 fVsafe(NULL),
120 fBuff1(NULL),
121 fBuff2(NULL),
122 fFullCheck(kFALSE),
123 fVal1(NULL),
124 fVal2(NULL),
125 fFlags(NULL),
126 fTimer(NULL),
127 fSelectedNode(NULL),
128 fNchecks(0),
129 fNmeshPoints(1000)
130{
131 fBuff1 = new TBuffer3D(TBuffer3DTypes::kGeneric,500,3*500,0,0,0,0);
132 fBuff2 = new TBuffer3D(TBuffer3DTypes::kGeneric,500,3*500,0,0,0,0);
133}
134
135////////////////////////////////////////////////////////////////////////////////
136/// Destructor
137
139{
140 if (fBuff1) delete fBuff1;
141 if (fBuff2) delete fBuff2;
142 if (fTimer) delete fTimer;
143}
144
145////////////////////////////////////////////////////////////////////////////////
146/// Print current operation progress.
147
148void TGeoChecker::OpProgress(const char *opname, Long64_t current, Long64_t size, TStopwatch *watch, Bool_t last, Bool_t refresh, const char *msg)
149{
150 static Long64_t icount = 0;
151 static TString oname;
152 static TString nname;
153 static Long64_t ocurrent = 0;
154 static Long64_t osize = 0;
155 static Int_t oseconds = 0;
156 static TStopwatch *owatch = 0;
157 static Bool_t oneoftwo = kFALSE;
158 static Int_t nrefresh = 0;
159 const char symbol[4] = {'=','\\','|','/'};
160 char progress[11] = " ";
161 Int_t ichar = icount%4;
162 TString message(msg);
163 message += " ";
164
165 if (!refresh) {
166 nrefresh = 0;
167 if (!size) return;
168 owatch = watch;
169 oname = opname;
170 ocurrent = TMath::Abs(current);
171 osize = TMath::Abs(size);
172 if (ocurrent > osize) ocurrent=osize;
173 } else {
174 nrefresh++;
175 if (!osize) return;
176 }
177 icount++;
178 Double_t time = 0.;
179 Int_t hours = 0;
180 Int_t minutes = 0;
181 Int_t seconds = 0;
182 if (owatch && !last) {
183 owatch->Stop();
184 time = owatch->RealTime();
185 hours = (Int_t)(time/3600.);
186 time -= 3600*hours;
187 minutes = (Int_t)(time/60.);
188 time -= 60*minutes;
189 seconds = (Int_t)time;
190 if (refresh) {
191 if (oseconds==seconds) {
192 owatch->Continue();
193 return;
194 }
195 oneoftwo = !oneoftwo;
196 }
197 oseconds = seconds;
198 }
199 if (refresh && oneoftwo) {
200 nname = oname;
201 if (fNchecks <= nrefresh) fNchecks = nrefresh+1;
202 Int_t pctdone = (Int_t)(100.*nrefresh/fNchecks);
203 oname = TString::Format(" == %3d%% ==", pctdone);
204 }
205 Double_t percent = 100.0*ocurrent/osize;
206 Int_t nchar = Int_t(percent/10);
207 if (nchar>10) nchar=10;
208 Int_t i;
209 for (i=0; i<nchar; i++) progress[i] = '=';
210 progress[nchar] = symbol[ichar];
211 for (i=nchar+1; i<10; i++) progress[i] = ' ';
212 progress[10] = '\0';
213 oname += " ";
214 oname.Remove(20);
215 if(size<10000) fprintf(stderr, "%s [%10s] %4lld ", oname.Data(), progress, ocurrent);
216 else if(size<100000) fprintf(stderr, "%s [%10s] %5lld ",oname.Data(), progress, ocurrent);
217 else fprintf(stderr, "%s [%10s] %7lld ",oname.Data(), progress, ocurrent);
218 if (time>0.) fprintf(stderr, "[%6.2f %%] TIME %.2d:%.2d:%.2d %s\r", percent, hours, minutes, seconds, message.Data());
219 else fprintf(stderr, "[%6.2f %%] %s\r", percent, message.Data());
220 if (refresh && oneoftwo) oname = nname;
221 if (owatch) owatch->Continue();
222 if (last) {
223 icount = 0;
224 owatch = 0;
225 ocurrent = 0;
226 osize = 0;
227 oseconds = 0;
228 oneoftwo = kFALSE;
229 nrefresh = 0;
230 fprintf(stderr, "\n");
231 }
232}
233
234////////////////////////////////////////////////////////////////////////////////
235/// Check pushes and pulls needed to cross the next boundary with respect to the
236/// position given by FindNextBoundary. If radius is not mentioned the full bounding
237/// box will be sampled.
238
240{
242 Info("CheckBoundaryErrors", "Top volume is %s",tvol->GetName());
243 const TGeoShape *shape = tvol->GetShape();
244 TGeoBBox *box = (TGeoBBox *)shape;
245 Double_t dl[3];
246 Double_t ori[3];
247 Double_t xyz[3];
248 Double_t nxyz[3];
249 Double_t dir[3];
250 Double_t relp;
251 Char_t path[1024];
252 Char_t cdir[10];
253
254 // Tree part
255 TFile *f=new TFile("geobugs.root","recreate");
256 TTree *bug=new TTree("bug","Geometrical problems");
257 bug->Branch("pos",xyz,"xyz[3]/D");
258 bug->Branch("dir",dir,"dir[3]/D");
259 bug->Branch("push",&relp,"push/D");
260 bug->Branch("path",&path,"path/C");
261 bug->Branch("cdir",&cdir,"cdir/C");
262
263 dl[0] = box->GetDX();
264 dl[1] = box->GetDY();
265 dl[2] = box->GetDZ();
266 ori[0] = (box->GetOrigin())[0];
267 ori[1] = (box->GetOrigin())[1];
268 ori[2] = (box->GetOrigin())[2];
269 if (radius>0)
270 dl[0] = dl[1] = dl[2] = radius;
271
273 TH1F *hnew = new TH1F("hnew","Precision pushing",30,-20.,10.);
274 TH1F *hold = new TH1F("hold","Precision pulling", 30,-20.,10.);
275 TH2F *hplotS = new TH2F("hplotS","Problematic points",100,-dl[0],dl[0],100,-dl[1],dl[1]);
276 gStyle->SetOptStat(111111);
277
278 TGeoNode *node = 0;
279 Long_t igen=0;
280 Long_t itry=0;
281 Long_t n100 = ntracks/100;
282 Double_t rad = TMath::Sqrt(dl[0]*dl[0]+dl[1]*dl[1]);
283 printf("Random box : %f, %f, %f, %f, %f, %f\n", ori[0], ori[1], ori[2], dl[0], dl[1], dl[2]);
284 printf("Start... %i points\n", ntracks);
285 if (!fTimer) fTimer = new TStopwatch();
286 fTimer->Reset();
287 fTimer->Start();
288 while (igen<ntracks) {
289 Double_t phi1 = TMath::TwoPi()*gRandom->Rndm();
290 Double_t r = rad*gRandom->Rndm();
291 xyz[0] = ori[0] + r*TMath::Cos(phi1);
292 xyz[1] = ori[1] + r*TMath::Sin(phi1);
293 Double_t z = (1.-2.*gRandom->Rndm());
294 xyz[2] = ori[2]+dl[2]*z*TMath::Abs(z);
295 ++itry;
297 node = fGeoManager->FindNode();
298 if (!node || node==fGeoManager->GetTopNode()) continue;
299 ++igen;
300 if (n100 && !(igen%n100))
301 OpProgress("Sampling progress:",igen, ntracks, fTimer);
302 Double_t cost = 1.-2.*gRandom->Rndm();
303 Double_t sint = TMath::Sqrt((1.+cost)*(1.-cost));
305 dir[0] = sint * TMath::Cos(phi);
306 dir[1] = sint * TMath::Sin(phi);
307 dir[2] = cost;
310 Double_t step = fGeoManager->GetStep();
311
312 relp = 1.e-21;
313 for(Int_t i=0; i<30; ++i) {
314 relp *=10.;
315 for(Int_t j=0; j<3; ++j) nxyz[j]=xyz[j]+step*(1.+relp)*dir[j];
316 if(!fGeoManager->IsSameLocation(nxyz[0],nxyz[1],nxyz[2])) {
317 hnew->Fill(i-20.);
318 if(i>15) {
319 const Double_t* norm = fGeoManager->FindNormal();
320 strncpy(path,fGeoManager->GetPath(),1024);
321 path[1023] = '\0';
322 Double_t dotp = norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2];
323 printf("Forward error i=%d p=%5.4f %5.4f %5.4f s=%5.4f dot=%5.4f path=%s\n",
324 i,xyz[0],xyz[1],xyz[2],step,dotp,path);
325 hplotS->Fill(xyz[0],xyz[1],(Double_t)i);
326 strncpy(cdir,"Forward",10);
327 bug->Fill();
328 }
329 break;
330 }
331 }
332
333 relp = -1.e-21;
334 for(Int_t i=0; i<30; ++i) {
335 relp *=10.;
336 for(Int_t j=0; j<3; ++j) nxyz[j]=xyz[j]+step*(1.+relp)*dir[j];
337 if(fGeoManager->IsSameLocation(nxyz[0],nxyz[1],nxyz[2])) {
338 hold->Fill(i-20.);
339 if(i>15) {
340 const Double_t* norm = fGeoManager->FindNormal();
341 strncpy(path,fGeoManager->GetPath(),1024);
342 path[1023] = '\0';
343 Double_t dotp = norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2];
344 printf("Backward error i=%d p=%5.4f %5.4f %5.4f s=%5.4f dot=%5.4f path=%s\n",
345 i,xyz[0],xyz[1],xyz[2],step,dotp,path);
346 strncpy(cdir,"Backward",10);
347 bug->Fill();
348 }
349 break;
350 }
351 }
352 }
353 fTimer->Stop();
354
355 if (itry) printf("CPU time/point = %5.2emus: Real time/point = %5.2emus\n",
356 1000000.*fTimer->CpuTime()/itry,1000000.*fTimer->RealTime()/itry);
357 bug->Write();
358 delete bug;
359 bug=0;
360 delete f;
361
363
364 if (itry) printf("Effic = %3.1f%%\n",(100.*igen)/itry);
365 TCanvas *c1 = new TCanvas("c1","Results",600,800);
366 c1->Divide(1,2);
367 c1->cd(1);
368 gPad->SetLogy();
369 hold->Draw();
370 c1->cd(2);
371 gPad->SetLogy();
372 hnew->Draw();
373 /*TCanvas *c3 = */new TCanvas("c3","Plot",600,600);
374 hplotS->Draw("cont0");
375}
376
377////////////////////////////////////////////////////////////////////////////////
378/// Check the boundary errors reference file created by CheckBoundaryErrors method.
379/// The shape for which the crossing failed is drawn with the starting point in red
380/// and the extrapolated point to boundary (+/- failing push/pull) in yellow.
381
383{
384 Double_t xyz[3];
385 Double_t nxyz[3];
386 Double_t dir[3];
387 Double_t lnext[3];
388 Double_t push;
389 Char_t path[1024];
390 Char_t cdir[10];
391 // Tree part
392 TFile *f=new TFile("geobugs.root","read");
393 TTree *bug=(TTree*)f->Get("bug");
394 bug->SetBranchAddress("pos",xyz);
395 bug->SetBranchAddress("dir",dir);
396 bug->SetBranchAddress("push",&push);
397 bug->SetBranchAddress("path",&path);
398 bug->SetBranchAddress("cdir",&cdir);
399
400 Int_t nentries = (Int_t)bug->GetEntries();
401 printf("nentries %d\n",nentries);
402 if (icheck<0) {
403 for (Int_t i=0;i<nentries;i++) {
404 bug->GetEntry(i);
405 printf("%-9s error push=%g p=%5.4f %5.4f %5.4f s=%5.4f dot=%5.4f path=%s\n",
406 cdir,push,xyz[0],xyz[1],xyz[2],1.,1.,path);
407 }
408 } else {
409 if (icheck>=nentries) return;
412 bug->GetEntry(icheck);
413 printf("%-9s error push=%g p=%5.4f %5.4f %5.4f s=%5.4f dot=%5.4f path=%s\n",
414 cdir,push,xyz[0],xyz[1],xyz[2],1.,1.,path);
419 Double_t step = fGeoManager->GetStep();
420 for (Int_t j=0; j<3; j++) nxyz[j]=xyz[j]+step*(1.+0.1*push)*dir[j];
421 Bool_t change = !fGeoManager->IsSameLocation(nxyz[0],nxyz[1],nxyz[2]);
422 printf("step=%g in: %s\n", step, fGeoManager->GetPath());
423 printf(" -> next = %s push=%g change=%d\n", next->GetName(),push, (UInt_t)change);
424 next->GetVolume()->InspectShape();
425 next->GetVolume()->DrawOnly();
426 if (next != fGeoManager->GetCurrentNode()) {
427 Int_t index1 = fGeoManager->GetCurrentVolume()->GetIndex(next);
428 if (index1>=0) fGeoManager->CdDown(index1);
429 }
430 TPolyMarker3D *pm = new TPolyMarker3D();
432 pm->SetNextPoint(lnext[0], lnext[1], lnext[2]);
433 pm->SetMarkerStyle(2);
434 pm->SetMarkerSize(0.2);
435 pm->SetMarkerColor(kRed);
436 pm->Draw("SAME");
437 TPolyMarker3D *pm1 = new TPolyMarker3D();
438 for (Int_t j=0; j<3; j++) nxyz[j]=xyz[j]+step*dir[j];
440 pm1->SetNextPoint(lnext[0], lnext[1], lnext[2]);
441 pm1->SetMarkerStyle(2);
442 pm1->SetMarkerSize(0.2);
444 pm1->Draw("SAME");
446 }
447 delete bug;
448 delete f;
449}
450
451////////////////////////////////////////////////////////////////////////////////
452/// Geometry checking. Optional overlap checkings (by sampling and by mesh). Optional
453/// boundary crossing check + timing per volume.
454///
455/// STAGE 1: extensive overlap checking by sampling per volume. Stdout need to be
456/// checked by user to get report, then TGeoVolume::CheckOverlaps(0.01, "s") can
457/// be called for the suspicious volumes.
458///
459/// STAGE 2: normal overlap checking using the shapes mesh - fills the list of
460/// overlaps.
461///
462/// STAGE 3: shooting NRAYS rays from VERTEX and counting the total number of
463/// crossings per volume (rays propagated from boundary to boundary until
464/// geometry exit). Timing computed and results stored in a histo.
465///
466/// STAGE 4: shooting 1 mil. random rays inside EACH volume and calling
467/// FindNextBoundary() + Safety() for each call. The timing is normalized by the
468/// number of crossings computed at stage 2 and presented as percentage.
469/// One can get a picture on which are the most "burned" volumes during
470/// transportation from geometry point of view. Another plot of the timing per
471/// volume vs. number of daughters is produced.
472///
473/// All histos are saved in the file statistics.root
474
475void TGeoChecker::CheckGeometryFull(Bool_t checkoverlaps, Bool_t checkcrossings, Int_t ntracks, const Double_t *vertex)
476{
478 if (!fTimer) fTimer = new TStopwatch();
479 Int_t i;
481 fFlags = new Bool_t[nuid];
482 memset(fFlags, 0, nuid*sizeof(Bool_t));
483 TGeoVolume *vol;
484 TCanvas *c = new TCanvas("overlaps", "Overlaps by sampling", 800,800);
485
486// STAGE 1: Overlap checking by sampling per volume
487 if (checkoverlaps) {
488 printf("====================================================================\n");
489 printf("STAGE 1: Overlap checking by sampling within 10 microns\n");
490 printf("====================================================================\n");
491 fGeoManager->CheckOverlaps(0.001, "s");
492
493 // STAGE 2: Global overlap/extrusion checking
494 printf("====================================================================\n");
495 printf("STAGE 2: Global overlap/extrusion checking within 10 microns\n");
496 printf("====================================================================\n");
498 }
499
500 if (!checkcrossings) {
501 delete [] fFlags;
502 fFlags = 0;
503 delete c;
504 return;
505 }
506
507 fVal1 = new Double_t[nuid];
508 fVal2 = new Double_t[nuid];
509 memset(fVal1, 0, nuid*sizeof(Double_t));
510 memset(fVal2, 0, nuid*sizeof(Double_t));
511 // STAGE 3: How many crossings per volume in a realistic case ?
512 // Ignore volumes with no daughters
513
514 // Generate rays from vertex in phi=[0,2*pi] theta=[0,pi]
515// Int_t ntracks = 1000000;
516 printf("====================================================================\n");
517 printf("STAGE 3: Propagating %i tracks starting from vertex\n and counting number of boundary crossings...\n", ntracks);
518 printf("====================================================================\n");
519 Int_t nbound = 0;
520 Double_t theta, phi;
521 Double_t point[3], dir[3];
522 memset(point, 0, 3*sizeof(Double_t));
523 if (vertex) memcpy(point, vertex, 3*sizeof(Double_t));
524
525 fTimer->Reset();
526 fTimer->Start();
527 for (i=0; i<ntracks; i++) {
528 phi = 2.*TMath::Pi()*gRandom->Rndm();
529 theta= TMath::ACos(1.-2.*gRandom->Rndm());
530 dir[0]=TMath::Sin(theta)*TMath::Cos(phi);
531 dir[1]=TMath::Sin(theta)*TMath::Sin(phi);
532 dir[2]=TMath::Cos(theta);
533 if ((i%100)==0) OpProgress("Transporting tracks",i, ntracks, fTimer);
534// if ((i%1000)==0) printf("... remaining tracks %i\n", ntracks-i);
535 nbound += PropagateInGeom(point,dir);
536 }
537 if (!nbound) {
538 printf("No boundary crossed\n");
539 return;
540 }
541 fTimer->Stop();
542 Double_t time1 = fTimer->CpuTime() *1.E6;
543 Double_t time2 = time1/ntracks;
544 Double_t time3 = time1/nbound;
545 OpProgress("Transporting tracks",ntracks, ntracks, fTimer, kTRUE);
546 printf("Time for crossing %i boundaries: %g [ms]\n", nbound, time1);
547 printf("Time per track for full geometry traversal: %g [ms], per crossing: %g [ms]\n", time2, time3);
548
549// STAGE 4: How much time per volume:
550
551 printf("====================================================================\n");
552 printf("STAGE 4: How much navigation time per volume per next+safety call\n");
553 printf("====================================================================\n");
555 TGeoNode*current;
556 TString path;
557 vol = fGeoManager->GetTopVolume();
558 memset(fFlags, 0, nuid*sizeof(Bool_t));
559 TStopwatch timer;
560 timer.Start();
561 i = 0;
562 char volname[30];
563 strncpy(volname, vol->GetName(),15);
564 volname[15] = '\0';
565 OpProgress(volname,i++, nuid, &timer);
566 Score(vol, 1, TimingPerVolume(vol));
567 while ((current=next())) {
568 vol = current->GetVolume();
569 Int_t uid = vol->GetNumber();
570 if (fFlags[uid]) continue;
571 fFlags[uid] = kTRUE;
572 next.GetPath(path);
573 fGeoManager->cd(path.Data());
574 strncpy(volname, vol->GetName(),15);
575 volname[15] = '\0';
576 OpProgress(volname,i++, nuid, &timer);
577 Score(vol,1,TimingPerVolume(vol));
578 }
579 OpProgress("STAGE 4 completed",i, nuid, &timer, kTRUE);
580 // Draw some histos
581 Double_t time_tot_pertrack = 0.;
582 TCanvas *c1 = new TCanvas("c2","ncrossings",10,10,900,500);
583 c1->SetGrid();
584 c1->SetTopMargin(0.15);
585 TFile *f = new TFile("statistics.root", "RECREATE");
586 TH1F *h = new TH1F("h","number of boundary crossings per volume",3,0,3);
587 h->SetStats(0);
588 h->SetFillColor(38);
589 h->SetCanExtend(TH1::kAllAxes);
590
591 memset(fFlags, 0, nuid*sizeof(Bool_t));
592 for (i=0; i<nuid; i++) {
593 vol = fGeoManager->GetVolume(i);
594 if (!vol->GetNdaughters()) continue;
595 time_tot_pertrack += fVal1[i]*fVal2[i];
596 h->Fill(vol->GetName(), (Int_t)fVal1[i]);
597 }
598 time_tot_pertrack /= ntracks;
599 h->LabelsDeflate();
600 h->LabelsOption(">","X");
601 h->Draw();
602
603
604 TCanvas *c2 = new TCanvas("c3","time spent per volume in navigation",10,10,900,500);
605 c2->SetGrid();
606 c2->SetTopMargin(0.15);
607 TH2F *h2 = new TH2F("h2", "time per FNB call vs. ndaughters", 100, 0,100,100,0,15);
608 h2->SetStats(0);
609 h2->SetMarkerStyle(2);
610 TH1F *h1 = new TH1F("h1","percent of time spent per volume",3,0,3);
611 h1->SetStats(0);
612 h1->SetFillColor(38);
614 for (i=0; i<nuid; i++) {
615 vol = fGeoManager->GetVolume(i);
616 if (!vol || !vol->GetNdaughters()) continue;
617 value = fVal1[i]*fVal2[i]/ntracks/time_tot_pertrack;
618 h1->Fill(vol->GetName(), value);
619 h2->Fill(vol->GetNdaughters(), fVal2[i]);
620 }
621 h1->LabelsDeflate();
622 h1->LabelsOption(">","X");
623 h1->Draw();
624 TCanvas *c3 = new TCanvas("c4","timing vs. ndaughters",10,10,900,500);
625 c3->SetGrid();
626 c3->SetTopMargin(0.15);
627 h2->Draw();
628 f->Write();
629 delete [] fFlags;
630 fFlags = 0;
631 delete [] fVal1;
632 fVal1 = 0;
633 delete [] fVal2;
634 fVal2 = 0;
635 delete fTimer;
636 fTimer = 0;
637 delete c;
638}
639
640////////////////////////////////////////////////////////////////////////////////
641/// Propagate from START along DIR from boundary to boundary until exiting
642/// geometry. Fill array of hits.
643
645{
646 fGeoManager->InitTrack(start, dir);
647 TGeoNode *current = 0;
648 Int_t nzero = 0;
649 Int_t nhits = 0;
650 while (!fGeoManager->IsOutside()) {
651 if (nzero>3) {
652 // Problems in trying to cross a boundary
653 printf("Error in trying to cross boundary of %s\n", current->GetName());
654 return nhits;
655 }
657 if (!current || fGeoManager->IsOutside()) return nhits;
658 Double_t step = fGeoManager->GetStep();
659 if (step<2.*TGeoShape::Tolerance()) {
660 nzero++;
661 continue;
662 }
663 else nzero = 0;
664 // Generate the hit
665 nhits++;
666 TGeoVolume *vol = current->GetVolume();
667 Score(vol,0,1.);
668 Int_t iup = 1;
669 TGeoNode *mother = fGeoManager->GetMother(iup++);
670 while (mother && mother->GetVolume()->IsAssembly()) {
671 Score(mother->GetVolume(), 0, 1.);
672 mother = fGeoManager->GetMother(iup++);
673 }
674 }
675 return nhits;
676}
677
678////////////////////////////////////////////////////////////////////////////////
679/// Score a hit for VOL
680
682{
683 Int_t uid = vol->GetNumber();
684 switch (ifield) {
685 case 0:
686 fVal1[uid] += value;
687 break;
688 case 1:
689 fVal2[uid] += value;
690 }
691}
692
693////////////////////////////////////////////////////////////////////////////////
694/// Set number of points to be generated on the shape outline when checking for overlaps.
695
697 fNmeshPoints = npoints;
698 if (npoints<1000) {
699 Error("SetNmeshPoints", "Cannot allow less than 1000 points for checking - set to 1000");
700 fNmeshPoints = 1000;
701 }
702}
703
704////////////////////////////////////////////////////////////////////////////////
705/// Compute timing per "FindNextBoundary" + "Safety" call. Volume must be
706/// in the current path.
707
709{
710 fTimer->Reset();
711 const TGeoShape *shape = vol->GetShape();
712 TGeoBBox *box = (TGeoBBox *)shape;
713 Double_t dx = box->GetDX();
714 Double_t dy = box->GetDY();
715 Double_t dz = box->GetDZ();
716 Double_t ox = (box->GetOrigin())[0];
717 Double_t oy = (box->GetOrigin())[1];
718 Double_t oz = (box->GetOrigin())[2];
719 Double_t point[3], dir[3], lpt[3], ldir[3];
720 Double_t pstep = 0.;
721 pstep = TMath::Max(pstep,dz);
722 Double_t theta, phi;
723 Int_t idaughter;
724 fTimer->Start();
725 for (Int_t i=0; i<1000000; i++) {
726 lpt[0] = ox-dx+2*dx*gRandom->Rndm();
727 lpt[1] = oy-dy+2*dy*gRandom->Rndm();
728 lpt[2] = oz-dz+2*dz*gRandom->Rndm();
730 fGeoManager->SetCurrentPoint(point[0],point[1],point[2]);
731 phi = 2*TMath::Pi()*gRandom->Rndm();
732 theta= TMath::ACos(1.-2.*gRandom->Rndm());
733 ldir[0]=TMath::Sin(theta)*TMath::Cos(phi);
734 ldir[1]=TMath::Sin(theta)*TMath::Sin(phi);
735 ldir[2]=TMath::Cos(theta);
738 fGeoManager->SetStep(pstep);
740 // dist = TGeoShape::Big();
741 if (!vol->IsAssembly()) {
742 Bool_t inside = vol->Contains(lpt);
743 if (!inside) {
744 // dist = vol->GetShape()->DistFromOutside(lpt,ldir,3,pstep);
745 // if (dist>=pstep) continue;
746 } else {
747 vol->GetShape()->DistFromInside(lpt,ldir,3,pstep);
748 }
749
750 if (!vol->GetNdaughters()) vol->GetShape()->Safety(lpt, inside);
751 }
752 if (vol->GetNdaughters()) {
754 fGeoManager->FindNextDaughterBoundary(point,dir,idaughter,kFALSE);
755 }
756 }
757 fTimer->Stop();
758 Double_t time_per_track = fTimer->CpuTime();
759 Int_t uid = vol->GetNumber();
760 Int_t ncrossings = (Int_t)fVal1[uid];
761 if (!vol->GetNdaughters())
762 printf("Time for volume %s (shape=%s): %g [ms] ndaughters=%d ncross=%d\n", vol->GetName(), vol->GetShape()->GetName(), time_per_track, vol->GetNdaughters(), ncrossings);
763 else
764 printf("Time for volume %s (assemb=%d): %g [ms] ndaughters=%d ncross=%d\n", vol->GetName(), vol->IsAssembly(), time_per_track, vol->GetNdaughters(), ncrossings);
765 return time_per_track;
766}
767
768////////////////////////////////////////////////////////////////////////////////
769/// Shoot nrays with random directions from starting point (startx, starty, startz)
770/// in the reference frame of this volume. Track each ray until exiting geometry, then
771/// shoot backwards from exiting point and compare boundary crossing points.
772
773void TGeoChecker::CheckGeometry(Int_t nrays, Double_t startx, Double_t starty, Double_t startz) const
774{
775 Int_t i, j;
776 Double_t start[3], end[3];
777 Double_t dir[3];
778 Double_t dummy[3];
779 Double_t eps = 0.;
780 Double_t *array1 = new Double_t[3*1000];
781 Double_t *array2 = new Double_t[3*1000];
782 TObjArray *pma = new TObjArray();
783 TPolyMarker3D *pm;
784 pm = new TPolyMarker3D();
785 pm->SetMarkerColor(2); // error > eps
786 pm->SetMarkerStyle(8);
787 pm->SetMarkerSize(0.4);
788 pma->AddAt(pm, 0);
789 pm = new TPolyMarker3D();
790 pm->SetMarkerColor(4); // point not found
791 pm->SetMarkerStyle(8);
792 pm->SetMarkerSize(0.4);
793 pma->AddAt(pm, 1);
794 pm = new TPolyMarker3D();
795 pm->SetMarkerColor(6); // extra point back
796 pm->SetMarkerStyle(8);
797 pm->SetMarkerSize(0.4);
798 pma->AddAt(pm, 2);
799 Int_t nelem1, nelem2;
800 Int_t dim1=1000, dim2=1000;
801 if ((startx==0) && (starty==0) && (startz==0)) eps=1E-3;
802 start[0] = startx+eps;
803 start[1] = starty+eps;
804 start[2] = startz+eps;
805 Int_t n10=nrays/10;
806 Double_t theta,phi;
807 Double_t dw, dwmin, dx, dy, dz;
808 Int_t ist1, ist2, ifound;
809 for (i=0; i<nrays; i++) {
810 if (n10) {
811 if ((i%n10) == 0) printf("%i percent\n", Int_t(100*i/nrays));
812 }
813 phi = 2*TMath::Pi()*gRandom->Rndm();
814 theta= TMath::ACos(1.-2.*gRandom->Rndm());
815 dir[0]=TMath::Sin(theta)*TMath::Cos(phi);
816 dir[1]=TMath::Sin(theta)*TMath::Sin(phi);
817 dir[2]=TMath::Cos(theta);
818 // shoot direct ray
819 nelem1=nelem2=0;
820// printf("DIRECT %i\n", i);
821 array1 = ShootRay(&start[0], dir[0], dir[1], dir[2], array1, nelem1, dim1);
822 if (!nelem1) continue;
823// for (j=0; j<nelem1; j++) printf("%i : %f %f %f\n", j, array1[3*j], array1[3*j+1], array1[3*j+2]);
824 memcpy(&end[0], &array1[3*(nelem1-1)], 3*sizeof(Double_t));
825 // shoot ray backwards
826// printf("BACK %i\n", i);
827 array2 = ShootRay(&end[0], -dir[0], -dir[1], -dir[2], array2, nelem2, dim2, &start[0]);
828 if (!nelem2) {
829 printf("#### NOTHING BACK ###########################\n");
830 for (j=0; j<nelem1; j++) {
831 pm = (TPolyMarker3D*)pma->At(0);
832 pm->SetNextPoint(array1[3*j], array1[3*j+1], array1[3*j+2]);
833 }
834 continue;
835 }
836// printf("BACKWARDS\n");
837 Int_t k=nelem2>>1;
838 for (j=0; j<k; j++) {
839 memcpy(&dummy[0], &array2[3*j], 3*sizeof(Double_t));
840 memcpy(&array2[3*j], &array2[3*(nelem2-1-j)], 3*sizeof(Double_t));
841 memcpy(&array2[3*(nelem2-1-j)], &dummy[0], 3*sizeof(Double_t));
842 }
843// for (j=0; j<nelem2; j++) printf("%i : %f ,%f ,%f \n", j, array2[3*j], array2[3*j+1], array2[3*j+2]);
844 if (nelem1!=nelem2) printf("### DIFFERENT SIZES : nelem1=%i nelem2=%i ##########\n", nelem1, nelem2);
845 ist1 = ist2 = 0;
846 // check first match
847
848 dx = array1[3*ist1]-array2[3*ist2];
849 dy = array1[3*ist1+1]-array2[3*ist2+1];
850 dz = array1[3*ist1+2]-array2[3*ist2+2];
851 dw = dx*dir[0]+dy*dir[1]+dz*dir[2];
852 fGeoManager->SetCurrentPoint(&array1[3*ist1]);
854// printf("%i : %s (%g, %g, %g)\n", ist1, fGeoManager->GetPath(),
855// array1[3*ist1], array1[3*ist1+1], array1[3*ist1+2]);
856 if (TMath::Abs(dw)<1E-4) {
857// printf(" matching %i (%g, %g, %g)\n", ist2, array2[3*ist2], array2[3*ist2+1], array2[3*ist2+2]);
858 ist2++;
859 } else {
860 printf("### NOT MATCHING %i f:(%f, %f, %f) b:(%f %f %f) DCLOSE=%f\n", ist2, array1[3*ist1], array1[3*ist1+1], array1[3*ist1+2], array2[3*ist2], array2[3*ist2+1], array2[3*ist2+2],dw);
861 pm = (TPolyMarker3D*)pma->At(0);
862 pm->SetNextPoint(array2[3*ist2], array2[3*ist2+1], array2[3*ist2+2]);
863 if (dw<0) {
864 // first boundary missed on way back
865 } else {
866 // first boundary different on way back
867 ist2++;
868 }
869 }
870
871 while ((ist1<nelem1-1) && (ist2<nelem2)) {
872 fGeoManager->SetCurrentPoint(&array1[3*ist1+3]);
874// printf("%i : %s (%g, %g, %g)\n", ist1+1, fGeoManager->GetPath(),
875// array1[3*ist1+3], array1[3*ist1+4], array1[3*ist1+5]);
876
877 dx = array1[3*ist1+3]-array1[3*ist1];
878 dy = array1[3*ist1+4]-array1[3*ist1+1];
879 dz = array1[3*ist1+5]-array1[3*ist1+2];
880 // distance to next point
881 dwmin = dx+dir[0]+dy*dir[1]+dz*dir[2];
882 while (ist2<nelem2) {
883 ifound = 0;
884 dx = array2[3*ist2]-array1[3*ist1];
885 dy = array2[3*ist2+1]-array1[3*ist1+1];
886 dz = array2[3*ist2+2]-array1[3*ist1+2];
887 dw = dx+dir[0]+dy*dir[1]+dz*dir[2];
888 if (TMath::Abs(dw-dwmin)<1E-4) {
889 ist1++;
890 ist2++;
891 break;
892 }
893 if (dw<dwmin) {
894 // point found on way back. Check if close enough to ist1+1
895 ifound++;
896 dw = dwmin-dw;
897 if (dw<1E-4) {
898 // point is matching ist1+1
899// printf(" matching %i (%g, %g, %g) DCLOSE=%g\n", ist2, array2[3*ist2], array2[3*ist2+1], array2[3*ist2+2], dw);
900 ist2++;
901 ist1++;
902 break;
903 } else {
904 // extra boundary found on way back
905 fGeoManager->SetCurrentPoint(&array2[3*ist2]);
907 pm = (TPolyMarker3D*)pma->At(2);
908 pm->SetNextPoint(array2[3*ist2], array2[3*ist2+1], array2[3*ist2+2]);
909 printf("### EXTRA BOUNDARY %i : %s found at DCLOSE=%f\n", ist2, fGeoManager->GetPath(), dw);
910 ist2++;
911 continue;
912 }
913 } else {
914 if (!ifound) {
915 // point ist1+1 not found on way back
916 fGeoManager->SetCurrentPoint(&array1[3*ist1+3]);
918 pm = (TPolyMarker3D*)pma->At(1);
919 pm->SetNextPoint(array2[3*ist1+3], array2[3*ist1+4], array2[3*ist1+5]);
920 printf("### BOUNDARY MISSED BACK #########################\n");
921 ist1++;
922 break;
923 } else {
924 ist1++;
925 break;
926 }
927 }
928 }
929 }
930 }
931 pm = (TPolyMarker3D*)pma->At(0);
932 pm->Draw("SAME");
933 pm = (TPolyMarker3D*)pma->At(1);
934 pm->Draw("SAME");
935 pm = (TPolyMarker3D*)pma->At(2);
936 pm->Draw("SAME");
937 if (gPad) {
938 gPad->Modified();
939 gPad->Update();
940 }
941 delete [] array1;
942 delete [] array2;
943 delete pma; // markers are drawn on the pad
944}
945
946////////////////////////////////////////////////////////////////////////////////
947/// Clean-up the mesh of pcon/pgon from useless points
948
950{
951 Int_t ipoint = 0;
952 Int_t j, k=0;
953 Double_t rsq;
954 for (Int_t i=0; i<numPoints; i++) {
955 j = 3*i;
956 rsq = points[j]*points[j]+points[j+1]*points[j+1];
957 if (rsq < 1.e-10) continue;
958 points[k] = points[j];
959 points[k+1] = points[j+1];
960 points[k+2] = points[j+2];
961 ipoint++;
962 k = 3*ipoint;
963 }
964 numPoints = ipoint;
965}
966
967////////////////////////////////////////////////////////////////////////////////
968/// Check if the 2 non-assembly volume candidates overlap/extrude. Returns overlap object.
969
971{
972 TGeoOverlap *nodeovlp = nullptr;
973 if (vol1->IsAssembly() || vol2->IsAssembly())
974 return nullptr;
975 if (dynamic_cast<TGeoTessellated *>(vol1->GetShape()) || dynamic_cast<TGeoTessellated *>(vol2->GetShape()))
976 return nullptr;
977 Int_t numPoints1 = fBuff1->NbPnts();
978 Int_t numSegs1 = fBuff1->NbSegs();
979 Int_t numPols1 = fBuff1->NbPols();
980 Int_t numPoints2 = fBuff2->NbPnts();
981 Int_t numSegs2 = fBuff2->NbSegs();
982 Int_t numPols2 = fBuff2->NbPols();
983 Int_t ip;
984 Bool_t extrude, isextrusion, isoverlapping;
985 Double_t *points1 = fBuff1->fPnts;
986 Double_t *points2 = fBuff2->fPnts;
987 Double_t local[3], local1[3];
988 Double_t point[3];
989 Double_t safety = TGeoShape::Big();
990 Double_t tolerance = TGeoShape::Tolerance();
991 TGeoShape *shape1 = vol1->GetShape();
992 TGeoShape *shape2 = vol2->GetShape();
993 OpProgress("refresh", 0,0,NULL,kFALSE,kTRUE);
994 shape1->GetMeshNumbers(numPoints1, numSegs1, numPols1);
995 if (fBuff1->fID != (TObject*)shape1) {
996 // Fill first buffer.
997 fBuff1->SetRawSizes(TMath::Max(numPoints1,fNmeshPoints), 3*TMath::Max(numPoints1,fNmeshPoints), 0, 0, 0, 0);
998 points1 = fBuff1->fPnts;
999 if (shape1->GetPointsOnSegments(fNmeshPoints, points1)) {
1000 numPoints1 = fNmeshPoints;
1001 } else {
1002 shape1->SetPoints(points1);
1003 }
1004// if (shape1->InheritsFrom(TGeoPcon::Class())) {
1005// CleanPoints(points1, numPoints1);
1006// fBuff1->SetRawSizes(numPoints1, 3*numPoints1,0, 0, 0, 0);
1007// }
1008 fBuff1->fID = shape1;
1009 }
1010 shape2->GetMeshNumbers(numPoints2, numSegs2, numPols2);
1011 if (fBuff2->fID != (TObject*)shape2) {
1012 // Fill second buffer.
1013 fBuff2->SetRawSizes(TMath::Max(numPoints2,fNmeshPoints), 3*TMath::Max(numPoints2,fNmeshPoints), 0, 0, 0,0);
1014 points2 = fBuff2->fPnts;
1015 if (shape2->GetPointsOnSegments(fNmeshPoints, points2)) {
1016 numPoints2 = fNmeshPoints;
1017 } else {
1018 shape2->SetPoints(points2);
1019 }
1020// if (shape2->InheritsFrom(TGeoPcon::Class())) {
1021// CleanPoints(points2, numPoints2);
1022// fBuff1->SetRawSizes(numPoints2, 3*numPoints2,0,0,0,0);
1023// }
1024 fBuff2->fID = shape2;
1025 }
1026
1027 if (!isovlp) {
1028 // Extrusion case. Test vol2 extrude vol1.
1029 isextrusion=kFALSE;
1030 // loop all points of the daughter
1031 for (ip=0; ip<numPoints2; ip++) {
1032 memcpy(local, &points2[3*ip], 3*sizeof(Double_t));
1033 if (TMath::Abs(local[0])<tolerance && TMath::Abs(local[1])<tolerance) continue;
1034 mat2->LocalToMaster(local, point);
1035 mat1->MasterToLocal(point, local);
1036 extrude = !shape1->Contains(local);
1037 if (extrude) {
1038 safety = shape1->Safety(local, kFALSE);
1039 if (safety<ovlp) extrude=kFALSE;
1040 }
1041 if (extrude) {
1042 if (!isextrusion) {
1043 isextrusion = kTRUE;
1044 nodeovlp = new TGeoOverlap(name, vol1, vol2, mat1,mat2,kFALSE,safety);
1045 nodeovlp->SetNextPoint(point[0],point[1],point[2]);
1046 fGeoManager->AddOverlap(nodeovlp);
1047 } else {
1048 if (safety>nodeovlp->GetOverlap()) nodeovlp->SetOverlap(safety);
1049 nodeovlp->SetNextPoint(point[0],point[1],point[2]);
1050 }
1051 }
1052 }
1053 // loop all points of the mother
1054 for (ip=0; ip<numPoints1; ip++) {
1055 memcpy(local, &points1[3*ip], 3*sizeof(Double_t));
1056 if (local[0]<1e-10 && local[1]<1e-10) continue;
1057 mat1->LocalToMaster(local, point);
1058 mat2->MasterToLocal(point, local1);
1059 extrude = shape2->Contains(local1);
1060 if (extrude) {
1061 // skip points on mother mesh that have no neighborhood outside mother
1062 safety = shape1->Safety(local,kTRUE);
1063 if (safety>1E-6) {
1064 extrude = kFALSE;
1065 } else {
1066 safety = shape2->Safety(local1,kTRUE);
1067 if (safety<ovlp) extrude=kFALSE;
1068 }
1069 }
1070 if (extrude) {
1071 if (!isextrusion) {
1072 isextrusion = kTRUE;
1073 nodeovlp = new TGeoOverlap(name, vol1,vol2,mat1,mat2,kFALSE,safety);
1074 nodeovlp->SetNextPoint(point[0],point[1],point[2]);
1075 fGeoManager->AddOverlap(nodeovlp);
1076 } else {
1077 if (safety>nodeovlp->GetOverlap()) nodeovlp->SetOverlap(safety);
1078 nodeovlp->SetNextPoint(point[0],point[1],point[2]);
1079 }
1080 }
1081 }
1082 return nodeovlp;
1083 }
1084 // Check overlap
1085 isoverlapping = kFALSE;
1086 // loop all points of first candidate
1087 for (ip=0; ip<numPoints1; ip++) {
1088 memcpy(local, &points1[3*ip], 3*sizeof(Double_t));
1089 if (local[0]<1e-10 && local[1]<1e-10) continue;
1090 mat1->LocalToMaster(local, point);
1091 mat2->MasterToLocal(point, local); // now point in local reference of second
1092 Bool_t overlap = shape2->Contains(local);
1093 if (overlap) {
1094 safety = shape2->Safety(local, kTRUE);
1095 if (safety<ovlp) overlap=kFALSE;
1096 }
1097 if (overlap) {
1098 if (!isoverlapping) {
1099 isoverlapping = kTRUE;
1100 nodeovlp = new TGeoOverlap(name,vol1,vol2,mat1,mat2,kTRUE,safety);
1101 nodeovlp->SetNextPoint(point[0],point[1],point[2]);
1102 fGeoManager->AddOverlap(nodeovlp);
1103 } else {
1104 if (safety>nodeovlp->GetOverlap()) nodeovlp->SetOverlap(safety);
1105 nodeovlp->SetNextPoint(point[0],point[1],point[2]);
1106 }
1107 }
1108 }
1109 // loop all points of second candidate
1110 for (ip=0; ip<numPoints2; ip++) {
1111 memcpy(local, &points2[3*ip], 3*sizeof(Double_t));
1112 if (local[0]<1e-10 && local[1]<1e-10) continue;
1113 mat2->LocalToMaster(local, point);
1114 mat1->MasterToLocal(point, local); // now point in local reference of first
1115 Bool_t overlap = shape1->Contains(local);
1116 if (overlap) {
1117 safety = shape1->Safety(local, kTRUE);
1118 if (safety<ovlp) overlap=kFALSE;
1119 }
1120 if (overlap) {
1121 if (!isoverlapping) {
1122 isoverlapping = kTRUE;
1123 nodeovlp = new TGeoOverlap(name,vol1,vol2,mat1,mat2,kTRUE,safety);
1124 nodeovlp->SetNextPoint(point[0],point[1],point[2]);
1125 fGeoManager->AddOverlap(nodeovlp);
1126 } else {
1127 if (safety>nodeovlp->GetOverlap()) nodeovlp->SetOverlap(safety);
1128 nodeovlp->SetNextPoint(point[0],point[1],point[2]);
1129 }
1130 }
1131 }
1132 return nodeovlp;
1133}
1134
1135////////////////////////////////////////////////////////////////////////////////
1136/// Check illegal overlaps for volume VOL within a limit OVLP by sampling npoints
1137/// inside the volume shape.
1138
1140{
1141 Int_t nd = vol->GetNdaughters();
1142 if (nd<2) return;
1143 TGeoVoxelFinder *voxels = vol->GetVoxels();
1144 if (!voxels) return;
1145 if (voxels->NeedRebuild()) {
1146 voxels->Voxelize();
1147 vol->FindOverlaps();
1148 }
1149 TGeoBBox *box = (TGeoBBox*)vol->GetShape();
1150 TGeoShape *shape;
1151 TGeoNode *node;
1152 Double_t dx = box->GetDX();
1153 Double_t dy = box->GetDY();
1154 Double_t dz = box->GetDZ();
1155 Double_t pt[3];
1156 Double_t local[3];
1157 Int_t *check_list = nullptr;
1158 Int_t ncheck = 0;
1159 const Double_t *orig = box->GetOrigin();
1160 Int_t ipoint = 0;
1161 Int_t itry = 0;
1162 Int_t iovlp = 0;
1163 Int_t id=0, id0=0, id1=0;
1164 Bool_t in, incrt;
1165 Double_t safe;
1166 TString name1 = "";
1167 TString name2 = "";
1168 TGeoOverlap **flags = nullptr;
1169 TGeoNode *node1, *node2;
1170 Int_t novlps = 0;
1171 TGeoHMatrix mat1, mat2;
1172// Int_t tid = TGeoManager::ThreadId();
1174 TGeoStateInfo &td = *nav->GetCache()->GetInfo();
1175 while (ipoint < npoints) {
1176 // Shoot randomly in the bounding box.
1177 pt[0] = orig[0] - dx + 2.*dx*gRandom->Rndm();
1178 pt[1] = orig[1] - dy + 2.*dy*gRandom->Rndm();
1179 pt[2] = orig[2] - dz + 2.*dz*gRandom->Rndm();
1180 if (!vol->Contains(pt)) {
1181 itry++;
1182 if (itry>10000 && !ipoint) {
1183 Error("CheckOverlapsBySampling", "No point inside volume!!! - aborting");
1184 break;
1185 }
1186 continue;
1187 }
1188 // Check if the point is inside one or more daughters
1189 in = kFALSE;
1190 ipoint++;
1191 check_list = voxels->GetCheckList(pt, ncheck, td);
1192 if (!check_list || ncheck<2) continue;
1193 for (id=0; id<ncheck; id++) {
1194 id0 = check_list[id];
1195 node = vol->GetNode(id0);
1196 // Ignore MANY's
1197 if (node->IsOverlapping()) continue;
1198 node->GetMatrix()->MasterToLocal(pt,local);
1199 shape = node->GetVolume()->GetShape();
1200 incrt = shape->Contains(local);
1201 if (!incrt) continue;
1202 if (!in) {
1203 in = kTRUE;
1204 id1 = id0;
1205 continue;
1206 }
1207 // The point is inside 2 or more daughters, check safety
1208 safe = shape->Safety(local, kTRUE);
1209// if (safe < ovlp) continue;
1210 // We really have found an overlap -> store the point in a container
1211 iovlp++;
1212 if (!novlps) {
1213 if (flags) delete [] flags; // should never happen
1214 flags = new TGeoOverlap*[nd*nd];
1215 memset(flags, 0, nd*nd*sizeof(TGeoOverlap*));
1216 }
1217 TGeoOverlap *nodeovlp = flags[nd*id1+id0];
1218 if (!nodeovlp) {
1219 novlps++;
1220 node1 = vol->GetNode(id1);
1221 name1 = node1->GetName();
1222 mat1 = node1->GetMatrix();
1224 while (cindex >= 0) {
1225 node1 = node1->GetVolume()->GetNode(cindex);
1226 name1 += TString::Format("/%s", node1->GetName());
1227 mat1.Multiply(node1->GetMatrix());
1228 cindex = node1->GetVolume()->GetCurrentNodeIndex();
1229 }
1230 node2 = vol->GetNode(id0);
1231 name2 = node2->GetName();
1232 mat2 = node2->GetMatrix();
1233 cindex = node2->GetVolume()->GetCurrentNodeIndex();
1234 while (cindex >= 0) {
1235 node2 = node2->GetVolume()->GetNode(cindex);
1236 name2 += TString::Format("/%s", node2->GetName());
1237 mat2.Multiply(node2->GetMatrix());
1238 cindex = node2->GetVolume()->GetCurrentNodeIndex();
1239 }
1240 nodeovlp = new TGeoOverlap(TString::Format("Volume %s: node %s overlapping %s",
1241 vol->GetName(), name1.Data(), name2.Data()), node1->GetVolume(),node2->GetVolume(),
1242 &mat1,&mat2, kTRUE, safe);
1243 flags[nd*id1+id0] = nodeovlp;
1244 fGeoManager->AddOverlap(nodeovlp);
1245 }
1246 // Max 100 points per marker
1247 if (nodeovlp->GetPolyMarker()->GetN()<100) nodeovlp->SetNextPoint(pt[0],pt[1],pt[2]);
1248 if (nodeovlp->GetOverlap()<safe) nodeovlp->SetOverlap(safe);
1249 }
1250 }
1251 nav->GetCache()->ReleaseInfo();
1252 if (flags) delete [] flags;
1253 if (!novlps) return;
1254 Double_t capacity = vol->GetShape()->Capacity();
1255 capacity *= Double_t(iovlp)/Double_t(npoints);
1256 Double_t err = 1./TMath::Sqrt(Double_t(iovlp));
1257 Info("CheckOverlapsBySampling", "#Found %d overlaps adding-up to %g +/- %g [cm3] for daughters of %s",
1258 novlps, capacity, err*capacity, vol->GetName());
1259}
1260
1261////////////////////////////////////////////////////////////////////////////////
1262/// Compute number of overlaps combinations to check per volume
1263
1265{
1266 if (vol->GetFinder()) return 0;
1267 UInt_t nd = vol->GetNdaughters();
1268 if (!nd) return 0;
1269 Bool_t is_assembly = vol->IsAssembly();
1270 TGeoIterator next1(vol);
1271 TGeoIterator next2(vol);
1272 Int_t nchecks = 0;
1273 TGeoNode *node;
1274 UInt_t id;
1275 if (!is_assembly) {
1276 while ((node=next1())) {
1277 if (node->IsOverlapping()) {
1278 next1.Skip();
1279 continue;
1280 }
1281 if (!node->GetVolume()->IsAssembly()) {
1282 nchecks++;
1283 next1.Skip();
1284 }
1285 }
1286 }
1287 // now check if the daughters overlap with each other
1288 if (nd<2) return nchecks;
1289 TGeoVoxelFinder *vox = vol->GetVoxels();
1290 if (!vox) return nchecks;
1291
1292
1293 TGeoNode *node1, *node01, *node02;
1294 Int_t novlp;
1295 Int_t *ovlps;
1296 Int_t ko;
1297 UInt_t io;
1298 for (id=0; id<nd; id++) {
1299 node01 = vol->GetNode(id);
1300 if (node01->IsOverlapping()) continue;
1301 vox->FindOverlaps(id);
1302 ovlps = node01->GetOverlaps(novlp);
1303 if (!ovlps) continue;
1304 for (ko=0; ko<novlp; ko++) { // loop all possible overlapping candidates
1305 io = ovlps[ko]; // (node1, shaped, matrix1, points, fBuff1)
1306 if (io<=id) continue;
1307 node02 = vol->GetNode(io);
1308 if (node02->IsOverlapping()) continue;
1309
1310 // We have to check node against node1, but they may be assemblies
1311 if (node01->GetVolume()->IsAssembly()) {
1312 next1.Reset(node01->GetVolume());
1313 while ((node=next1())) {
1314 if (!node->GetVolume()->IsAssembly()) {
1315 if (node02->GetVolume()->IsAssembly()) {
1316 next2.Reset(node02->GetVolume());
1317 while ((node1=next2())) {
1318 if (!node1->GetVolume()->IsAssembly()) {
1319 nchecks++;
1320 next2.Skip();
1321 }
1322 }
1323 } else {
1324 nchecks++;
1325 }
1326 next1.Skip();
1327 }
1328 }
1329 } else {
1330 // node not assembly
1331 if (node02->GetVolume()->IsAssembly()) {
1332 next2.Reset(node02->GetVolume());
1333 while ((node1=next2())) {
1334 if (!node1->GetVolume()->IsAssembly()) {
1335 nchecks++;
1336 next2.Skip();
1337 }
1338 }
1339 } else {
1340 // node1 also not assembly
1341 nchecks++;
1342 }
1343 }
1344 }
1345 node01->SetOverlaps(0,0);
1346 }
1347 return nchecks;
1348}
1349
1350////////////////////////////////////////////////////////////////////////////////
1351/// Check illegal overlaps for volume VOL within a limit OVLP.
1352
1354{
1355 if (vol->GetFinder()) return;
1356 UInt_t nd = vol->GetNdaughters();
1357 if (!nd) return;
1360 Bool_t sampling = kFALSE;
1361 TString opt(option);
1362 opt.ToLower();
1363 if (opt.Contains("s")) sampling = kTRUE;
1364 if (opt.Contains("f")) fFullCheck = kTRUE;
1365 else fFullCheck = kFALSE;
1366 if (sampling) {
1367 opt = opt.Strip(TString::kLeading, 's');
1368 Int_t npoints = atoi(opt.Data());
1369 if (!npoints) npoints = 1000000;
1370 CheckOverlapsBySampling((TGeoVolume*)vol, ovlp, npoints);
1371 return;
1372 }
1373 Bool_t is_assembly = vol->IsAssembly();
1374 TGeoIterator next1((TGeoVolume*)vol);
1375 TGeoIterator next2((TGeoVolume*)vol);
1376 TString path;
1377 // first, test if daughters extrude their container
1378 TGeoNode * node, *nodecheck;
1379 TGeoChecker *checker = (TGeoChecker*)this;
1380
1381// TGeoOverlap *nodeovlp = 0;
1382 UInt_t id;
1383 Int_t level;
1384// Check extrusion only for daughters of a non-assembly volume
1385 if (!is_assembly) {
1386 while ((node=next1())) {
1387 if (node->IsOverlapping()) {
1388 next1.Skip();
1389 continue;
1390 }
1391 if (!node->GetVolume()->IsAssembly()) {
1392 if (fSelectedNode) {
1393 // We have to check only overlaps of the selected node (or real daughters if an assembly)
1394 if ((fSelectedNode != node) && (!fSelectedNode->GetVolume()->IsAssembly())) {
1395 next1.Skip();
1396 continue;
1397 }
1398 if (node != fSelectedNode) {
1399 level = next1.GetLevel();
1400 while ((nodecheck=next1.GetNode(level--))) {
1401 if (nodecheck == fSelectedNode) break;
1402 }
1403 if (!nodecheck) {
1404 next1.Skip();
1405 continue;
1406 }
1407 }
1408 }
1409 next1.GetPath(path);
1410 checker->MakeCheckOverlap(TString::Format("%s extruded by: %s", vol->GetName(),path.Data()),
1411 (TGeoVolume*)vol,node->GetVolume(),gGeoIdentity,(TGeoMatrix*)next1.GetCurrentMatrix(),kFALSE,ovlp);
1412 next1.Skip();
1413 }
1414 }
1415 }
1416
1417 // now check if the daughters overlap with each other
1418 if (nd<2) return;
1419 TGeoVoxelFinder *vox = vol->GetVoxels();
1420 if (!vox) {
1421 Warning("CheckOverlaps", "Volume %s with %i daughters but not voxelized", vol->GetName(),nd);
1422 return;
1423 }
1424 if (vox->NeedRebuild()) {
1425 vox->Voxelize();
1426 vol->FindOverlaps();
1427 }
1428 TGeoNode *node1, *node01, *node02;
1429 TGeoHMatrix hmat1, hmat2;
1430 TString path1;
1431 Int_t novlp;
1432 Int_t *ovlps;
1433 Int_t ko;
1434 UInt_t io;
1435 for (id=0; id<nd; id++) {
1436 node01 = vol->GetNode(id);
1437 if (node01->IsOverlapping()) continue;
1438 vox->FindOverlaps(id);
1439 ovlps = node01->GetOverlaps(novlp);
1440 if (!ovlps) continue;
1441 next1.SetTopName(node01->GetName());
1442 path = node01->GetName();
1443 for (ko=0; ko<novlp; ko++) { // loop all possible overlapping candidates
1444 io = ovlps[ko]; // (node1, shaped, matrix1, points, fBuff1)
1445 if (io<=id) continue;
1446 node02 = vol->GetNode(io);
1447 if (node02->IsOverlapping()) continue;
1448 // Try to fasten-up things...
1449// if (!TGeoBBox::AreOverlapping((TGeoBBox*)node01->GetVolume()->GetShape(), node01->GetMatrix(),
1450// (TGeoBBox*)node02->GetVolume()->GetShape(), node02->GetMatrix())) continue;
1451 next2.SetTopName(node02->GetName());
1452 path1 = node02->GetName();
1453
1454 // We have to check node against node1, but they may be assemblies
1455 if (node01->GetVolume()->IsAssembly()) {
1456 next1.Reset(node01->GetVolume());
1457 while ((node=next1())) {
1458 if (!node->GetVolume()->IsAssembly()) {
1459 next1.GetPath(path);
1460 hmat1 = node01->GetMatrix();
1461 hmat1 *= *next1.GetCurrentMatrix();
1462 if (node02->GetVolume()->IsAssembly()) {
1463 next2.Reset(node02->GetVolume());
1464 while ((node1=next2())) {
1465 if (!node1->GetVolume()->IsAssembly()) {
1466 if (fSelectedNode) {
1467 // We have to check only overlaps of the selected node (or real daughters if an assembly)
1468 if ((fSelectedNode != node) && (fSelectedNode != node1) && (!fSelectedNode->GetVolume()->IsAssembly())) {
1469 next2.Skip();
1470 continue;
1471 }
1472 if ((fSelectedNode != node1) && (fSelectedNode != node)) {
1473 level = next2.GetLevel();
1474 while ((nodecheck=next2.GetNode(level--))) {
1475 if (nodecheck == fSelectedNode) break;
1476 }
1477 if (node02 == fSelectedNode) nodecheck = node02;
1478 if (!nodecheck) {
1479 level = next1.GetLevel();
1480 while ((nodecheck=next1.GetNode(level--))) {
1481 if (nodecheck == fSelectedNode) break;
1482 }
1483 }
1484 if (node01 == fSelectedNode) nodecheck = node01;
1485 if (!nodecheck) {
1486 next2.Skip();
1487 continue;
1488 }
1489 }
1490 }
1491 next2.GetPath(path1);
1492 hmat2 = node02->GetMatrix();
1493 hmat2 *= *next2.GetCurrentMatrix();
1494 checker->MakeCheckOverlap(TString::Format("%s/%s overlapping %s/%s", vol->GetName(),path.Data(),vol->GetName(),path1.Data()),
1495 node->GetVolume(),node1->GetVolume(),&hmat1,&hmat2,kTRUE,ovlp);
1496 next2.Skip();
1497 }
1498 }
1499 } else {
1500 if (fSelectedNode) {
1501 // We have to check only overlaps of the selected node (or real daughters if an assembly)
1502 if ((fSelectedNode != node) && (fSelectedNode != node02) && (!fSelectedNode->GetVolume()->IsAssembly())) {
1503 next1.Skip();
1504 continue;
1505 }
1506 if ((fSelectedNode != node) && (fSelectedNode != node02)) {
1507 level = next1.GetLevel();
1508 while ((nodecheck=next1.GetNode(level--))) {
1509 if (nodecheck == fSelectedNode) break;
1510 }
1511 if (node01 == fSelectedNode) nodecheck = node01;
1512 if (!nodecheck) {
1513 next1.Skip();
1514 continue;
1515 }
1516 }
1517 }
1518 checker->MakeCheckOverlap(TString::Format("%s/%s overlapping %s/%s", vol->GetName(),path.Data(),vol->GetName(),path1.Data()),
1519 node->GetVolume(),node02->GetVolume(),&hmat1,node02->GetMatrix(),kTRUE,ovlp);
1520 }
1521 next1.Skip();
1522 }
1523 }
1524 } else {
1525 // node not assembly
1526 if (node02->GetVolume()->IsAssembly()) {
1527 next2.Reset(node02->GetVolume());
1528 while ((node1=next2())) {
1529 if (!node1->GetVolume()->IsAssembly()) {
1530 if (fSelectedNode) {
1531 // We have to check only overlaps of the selected node (or real daughters if an assembly)
1532 if ((fSelectedNode != node1) && (fSelectedNode != node01) && (!fSelectedNode->GetVolume()->IsAssembly())) {
1533 next2.Skip();
1534 continue;
1535 }
1536 if ((fSelectedNode != node1) && (fSelectedNode != node01)) {
1537 level = next2.GetLevel();
1538 while ((nodecheck=next2.GetNode(level--))) {
1539 if (nodecheck == fSelectedNode) break;
1540 }
1541 if (node02 == fSelectedNode) nodecheck = node02;
1542 if (!nodecheck) {
1543 next2.Skip();
1544 continue;
1545 }
1546 }
1547 }
1548 next2.GetPath(path1);
1549 hmat2 = node02->GetMatrix();
1550 hmat2 *= *next2.GetCurrentMatrix();
1551 checker->MakeCheckOverlap(TString::Format("%s/%s overlapping %s/%s", vol->GetName(),path.Data(),vol->GetName(),path1.Data()),
1552 node01->GetVolume(),node1->GetVolume(),node01->GetMatrix(),&hmat2,kTRUE,ovlp);
1553 next2.Skip();
1554 }
1555 }
1556 } else {
1557 // node1 also not assembly
1558 if (fSelectedNode && (fSelectedNode != node01) && (fSelectedNode != node02)) continue;
1559 checker->MakeCheckOverlap(TString::Format("%s/%s overlapping %s/%s", vol->GetName(),path.Data(),vol->GetName(),path1.Data()),
1560 node01->GetVolume(),node02->GetVolume(),node01->GetMatrix(),node02->GetMatrix(),kTRUE,ovlp);
1561 }
1562 }
1563 }
1564 node01->SetOverlaps(0,0);
1565 }
1566}
1567
1568////////////////////////////////////////////////////////////////////////////////
1569/// Print the current list of overlaps held by the manager class.
1570
1572{
1574 TGeoOverlap *ov;
1575 printf("=== Overlaps for %s ===\n", fGeoManager->GetName());
1576 while ((ov=(TGeoOverlap*)next())) ov->PrintInfo();
1577}
1578
1579////////////////////////////////////////////////////////////////////////////////
1580/// Draw point (x,y,z) over the picture of the daughters of the volume containing this point.
1581/// Generates a report regarding the path to the node containing this point and the distance to
1582/// the closest boundary.
1583
1585{
1586 Double_t point[3];
1587 Double_t local[3];
1588 point[0] = x;
1589 point[1] = y;
1590 point[2] = z;
1592 if (fVsafe) {
1593 TGeoNode *old = fVsafe->GetNode("SAFETY_1");
1594 if (old) fVsafe->GetNodes()->RemoveAt(vol->GetNdaughters()-1);
1595 }
1596// if (vol != fGeoManager->GetMasterVolume()) fGeoManager->RestoreMasterVolume();
1597 TGeoNode *node = fGeoManager->FindNode(point[0], point[1], point[2]);
1598 fGeoManager->MasterToLocal(point, local);
1599 // get current node
1600 printf("=== Check current point : (%g, %g, %g) ===\n", point[0], point[1], point[2]);
1601 printf(" - path : %s\n", fGeoManager->GetPath());
1602 // get corresponding volume
1603 if (node) vol = node->GetVolume();
1604 // compute safety distance (distance to boundary ignored)
1605 Double_t close = fGeoManager->Safety();
1606 printf("Safety radius : %f\n", close);
1607 if (close>1E-4) {
1608 TGeoVolume *sph = fGeoManager->MakeSphere("SAFETY", vol->GetMedium(), 0, close, 0,180,0,360);
1609 sph->SetLineColor(2);
1610 sph->SetLineStyle(3);
1611 vol->AddNode(sph,1,new TGeoTranslation(local[0], local[1], local[2]));
1612 fVsafe = vol;
1613 }
1614 TPolyMarker3D *pm = new TPolyMarker3D();
1615 pm->SetMarkerColor(2);
1616 pm->SetMarkerStyle(8);
1617 pm->SetMarkerSize(0.5);
1618 pm->SetNextPoint(local[0], local[1], local[2]);
1619 if (vol->GetNdaughters()<2) fGeoManager->SetTopVisible();
1622 if (!vol->IsVisible()) vol->SetVisibility(kTRUE);
1623 vol->Draw();
1624 pm->Draw("SAME");
1625 gPad->Modified();
1626 gPad->Update();
1627}
1628
1629////////////////////////////////////////////////////////////////////////////////
1630/// Test for shape navigation methods. Summary for test numbers:
1631/// - 1: DistFromInside/Outside. Sample points inside the shape. Generate
1632/// directions randomly in cos(theta). Compute DistFromInside and move the
1633/// point with bigger distance. Compute DistFromOutside back from new point.
1634/// Plot d-(d1+d2)
1635/// - 2: Safety test. Sample points inside the bounding and compute safety. Generate
1636/// directions randomly in cos(theta) and compute distance to boundary. Check if
1637/// distance to boundary is bigger than safety
1638
1640{
1641 switch (testNo) {
1642 case 1:
1643 ShapeDistances(shape, nsamples, option);
1644 break;
1645 case 2:
1646 ShapeSafety(shape, nsamples, option);
1647 break;
1648 case 3:
1649 ShapeNormal(shape, nsamples, option);
1650 break;
1651 default:
1652 Error("CheckShape", "Test number %d not existent", testNo);
1653 }
1654}
1655
1656////////////////////////////////////////////////////////////////////////////////
1657/// Test TGeoShape::DistFromInside/Outside. Sample points inside the shape. Generate
1658/// directions randomly in cos(theta). Compute d1 = DistFromInside and move the
1659/// point on the boundary. Compute DistFromOutside and propagate with d2 making sure that
1660/// the shape is not re-entered. Swap direction and call DistFromOutside that
1661/// should fall back on the same point on the boundary (at d2). Propagate back on boundary
1662/// then compute DistFromInside that should be bigger than d1.
1663/// Plot d-(d1+d2)
1664
1666{
1667 Double_t dx = ((TGeoBBox*)shape)->GetDX();
1668 Double_t dy = ((TGeoBBox*)shape)->GetDY();
1669 Double_t dz = ((TGeoBBox*)shape)->GetDZ();
1670 Double_t dmax = 2.*TMath::Sqrt(dx*dx+dy*dy+dz*dz);
1671 Double_t d1, d2, dmove, dnext;
1672 Int_t itot = 0;
1673 // Number of tracks shot for every point inside the shape
1674 const Int_t kNtracks = 1000;
1675 Int_t n10 = nsamples/10;
1676 Int_t i,j;
1677 Double_t point[3], pnew[3];
1678 Double_t dir[3], dnew[3];
1679 Double_t theta, phi, delta;
1680 TPolyMarker3D *pmfrominside = 0;
1681 TPolyMarker3D *pmfromoutside = 0;
1682 TH1D *hist = new TH1D("hTest1", "Residual distance from inside/outside",200,-20, 0);
1683 hist->GetXaxis()->SetTitle("delta[cm] - first bin=overflow");
1684 hist->GetYaxis()->SetTitle("count");
1686
1687 if (!fTimer) fTimer = new TStopwatch();
1688 fTimer->Reset();
1689 fTimer->Start();
1690 while (itot<nsamples) {
1691 Bool_t inside = kFALSE;
1692 while (!inside) {
1693 point[0] = gRandom->Uniform(-dx,dx);
1694 point[1] = gRandom->Uniform(-dy,dy);
1695 point[2] = gRandom->Uniform(-dz,dz);
1696 inside = shape->Contains(point);
1697 }
1698 itot++;
1699 if (n10) {
1700 if ((itot%n10) == 0) printf("%i percent\n", Int_t(100*itot/nsamples));
1701 }
1702 for (i=0; i<kNtracks; i++) {
1703 phi = 2*TMath::Pi()*gRandom->Rndm();
1704 theta= TMath::ACos(1.-2.*gRandom->Rndm());
1705 dir[0]=TMath::Sin(theta)*TMath::Cos(phi);
1706 dir[1]=TMath::Sin(theta)*TMath::Sin(phi);
1707 dir[2]=TMath::Cos(theta);
1708 dmove = dmax;
1709 // We have track direction, compute distance from inside
1710 d1 = shape->DistFromInside(point,dir,3);
1711 if (d1>dmove || d1<TGeoShape::Tolerance()) {
1712 // Bad distance or bbox size, to debug
1713 new TCanvas("shape01", Form("Shape %s (%s)",shape->GetName(),shape->ClassName()), 1000, 800);
1714 shape->Draw();
1715 printf("DistFromInside: (%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) %f/%f(max)\n",
1716 point[0],point[1],point[2],dir[0],dir[1],dir[2], d1,dmove);
1717 pmfrominside = new TPolyMarker3D(2);
1718 pmfrominside->SetMarkerColor(kRed);
1719 pmfrominside->SetMarkerStyle(24);
1720 pmfrominside->SetMarkerSize(0.4);
1721 pmfrominside->SetNextPoint(point[0],point[1],point[2]);
1722 for (j=0; j<3; j++) pnew[j] = point[j] + d1*dir[j];
1723 pmfrominside->SetNextPoint(pnew[0],pnew[1],pnew[2]);
1724 pmfrominside->Draw();
1725 return;
1726 }
1727 // Propagate BEFORE the boundary and make sure that DistFromOutside
1728 // does not return 0 (!!!)
1729 // Check if there is a second crossing
1730 for (j=0; j<3; j++) pnew[j] = point[j] + (d1-TGeoShape::Tolerance())*dir[j];
1731 dnext = shape->DistFromOutside(pnew,dir,3);
1732 if (d1+dnext<dmax) dmove = d1+0.5*dnext;
1733 // Move point and swap direction
1734 for (j=0; j<3; j++) {
1735 pnew[j] = point[j] + dmove*dir[j];
1736 dnew[j] = -dir[j];
1737 }
1738 // Compute now distance from outside
1739 d2 = shape->DistFromOutside(pnew,dnew,3);
1740 delta = dmove-d1-d2;
1741 if (TMath::Abs(delta)>1E-6 || dnext<2.*TGeoShape::Tolerance()) {
1742 // Error->debug this
1743 new TCanvas("shape01", Form("Shape %s (%s)",shape->GetName(),shape->ClassName()), 1000, 800);
1744 shape->Draw();
1745 printf("Error: (%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) d1=%f d2=%f dmove=%f\n",
1746 point[0],point[1],point[2],dir[0],dir[1],dir[2], d1,d2,dmove);
1747 if (dnext<2.*TGeoShape::Tolerance()) {
1748 printf(" (*)DistFromOutside(%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) dnext = %f\n",
1749 point[0]+(d1-TGeoShape::Tolerance())*dir[0],
1750 point[1]+(d1-TGeoShape::Tolerance())*dir[1],
1751 point[2]+(d1-TGeoShape::Tolerance())*dir[2], dir[0],dir[1],dir[2],dnext);
1752 } else {
1753 printf(" DistFromOutside(%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) dnext = %f\n",
1754 point[0]+d1*dir[0],point[1]+d1*dir[1], point[2]+d1*dir[2], dir[0],dir[1],dir[2],dnext);
1755 }
1756 printf(" DistFromOutside(%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) = %f\n",
1757 pnew[0],pnew[1],pnew[2],dnew[0],dnew[1],dnew[2], d2);
1758 pmfrominside = new TPolyMarker3D(2);
1759 pmfrominside->SetMarkerStyle(24);
1760 pmfrominside->SetMarkerSize(0.4);
1761 pmfrominside->SetMarkerColor(kRed);
1762 pmfrominside->SetNextPoint(point[0],point[1],point[2]);
1763 for (j=0; j<3; j++) point[j] += d1*dir[j];
1764 pmfrominside->SetNextPoint(point[0],point[1],point[2]);
1765 pmfrominside->Draw();
1766 pmfromoutside = new TPolyMarker3D(2);
1767 pmfromoutside->SetMarkerStyle(20);
1768 pmfromoutside->SetMarkerStyle(7);
1769 pmfromoutside->SetMarkerSize(0.3);
1770 pmfromoutside->SetMarkerColor(kBlue);
1771 pmfromoutside->SetNextPoint(pnew[0],pnew[1],pnew[2]);
1772 for (j=0; j<3; j++) pnew[j] += d2*dnew[j];
1773 if (d2<1E10) pmfromoutside->SetNextPoint(pnew[0],pnew[1],pnew[2]);
1774 pmfromoutside->Draw();
1775 return;
1776 }
1777 // Compute distance from inside which should be bigger than d1
1778 for (j=0; j<3; j++) pnew[j] += (d2-TGeoShape::Tolerance())*dnew[j];
1779 dnext = shape->DistFromInside(pnew,dnew,3);
1780 if (dnext<d1-TGeoShape::Tolerance() || dnext>dmax) {
1781 new TCanvas("shape01", Form("Shape %s (%s)",shape->GetName(),shape->ClassName()), 1000, 800);
1782 shape->Draw();
1783 printf("Error DistFromInside(%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) d1=%f d1p=%f\n",
1784 pnew[0],pnew[1],pnew[2],dnew[0],dnew[1],dnew[2],d1,dnext);
1785 pmfrominside = new TPolyMarker3D(2);
1786 pmfrominside->SetMarkerStyle(24);
1787 pmfrominside->SetMarkerSize(0.4);
1788 pmfrominside->SetMarkerColor(kRed);
1789 pmfrominside->SetNextPoint(point[0],point[1],point[2]);
1790 for (j=0; j<3; j++) point[j] += d1*dir[j];
1791 pmfrominside->SetNextPoint(point[0],point[1],point[2]);
1792 pmfrominside->Draw();
1793 pmfromoutside = new TPolyMarker3D(2);
1794 pmfromoutside->SetMarkerStyle(20);
1795 pmfromoutside->SetMarkerStyle(7);
1796 pmfromoutside->SetMarkerSize(0.3);
1797 pmfromoutside->SetMarkerColor(kBlue);
1798 pmfromoutside->SetNextPoint(pnew[0],pnew[1],pnew[2]);
1799 for (j=0; j<3; j++) pnew[j] += dnext*dnew[j];
1800 if (d2<1E10) pmfromoutside->SetNextPoint(pnew[0],pnew[1],pnew[2]);
1801 pmfromoutside->Draw();
1802 return;
1803 }
1804 if (TMath::Abs(delta) < 1E-20) delta = 1E-30;
1805 hist->Fill(TMath::Max(TMath::Log(TMath::Abs(delta)),-20.));
1806 }
1807 }
1808 fTimer->Stop();
1809 fTimer->Print();
1810 new TCanvas("Test01", "Residuals DistFromInside/Outside", 800, 600);
1811 hist->Draw();
1812}
1813
1814////////////////////////////////////////////////////////////////////////////////
1815/// Check of validity of safe distance for a given shape.
1816/// Sample points inside the 2x bounding box and compute safety. Generate
1817/// directions randomly in cos(theta) and compute distance to boundary. Check if
1818/// distance to boundary is bigger than safety.
1819
1821{
1822 Double_t dx = ((TGeoBBox*)shape)->GetDX();
1823 Double_t dy = ((TGeoBBox*)shape)->GetDY();
1824 Double_t dz = ((TGeoBBox*)shape)->GetDZ();
1825 // Number of tracks shot for every point inside the shape
1826 const Int_t kNtracks = 1000;
1827 Int_t n10 = nsamples/10;
1828 Int_t i;
1829 Double_t dist;
1830 Double_t point[3];
1831 Double_t dir[3];
1832 Double_t theta, phi;
1833 TPolyMarker3D *pm1 = 0;
1834 TPolyMarker3D *pm2 = 0;
1835 if (!fTimer) fTimer = new TStopwatch();
1836 fTimer->Reset();
1837 fTimer->Start();
1838 Int_t itot = 0;
1839 while (itot<nsamples) {
1840 Bool_t inside = kFALSE;
1841 point[0] = gRandom->Uniform(-2*dx,2*dx);
1842 point[1] = gRandom->Uniform(-2*dy,2*dy);
1843 point[2] = gRandom->Uniform(-2*dz,2*dz);
1844 inside = shape->Contains(point);
1845 Double_t safe = shape->Safety(point, inside);
1846 itot++;
1847 if (n10) {
1848 if ((itot%n10) == 0) printf("%i percent\n", Int_t(100*itot/nsamples));
1849 }
1850 for (i=0; i<kNtracks; i++) {
1851 phi = 2*TMath::Pi()*gRandom->Rndm();
1852 theta= TMath::ACos(1.-2.*gRandom->Rndm());
1853 dir[0]=TMath::Sin(theta)*TMath::Cos(phi);
1854 dir[1]=TMath::Sin(theta)*TMath::Sin(phi);
1855 dir[2]=TMath::Cos(theta);
1856 if (inside) dist = shape->DistFromInside(point,dir,3);
1857 else dist = shape->DistFromOutside(point,dir,3);
1858 if (dist<safe) {
1859 printf("Error safety (%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) safe=%f dist=%f\n",
1860 point[0],point[1],point[2], dir[0], dir[1], dir[2], safe, dist);
1861 new TCanvas("shape02", Form("Shape %s (%s)",shape->GetName(),shape->ClassName()), 1000, 800);
1862 shape->Draw();
1863 pm1 = new TPolyMarker3D(2);
1864 pm1->SetMarkerStyle(24);
1865 pm1->SetMarkerSize(0.4);
1866 pm1->SetMarkerColor(kRed);
1867 pm1->SetNextPoint(point[0],point[1],point[2]);
1868 pm1->SetNextPoint(point[0]+safe*dir[0],point[1]+safe*dir[1],point[2]+safe*dir[2]);
1869 pm1->Draw();
1870 pm2 = new TPolyMarker3D(1);
1871 pm2->SetMarkerStyle(7);
1872 pm2->SetMarkerSize(0.3);
1873 pm2->SetMarkerColor(kBlue);
1874 pm2->SetNextPoint(point[0]+dist*dir[0],point[1]+dist*dir[1],point[2]+dist*dir[2]);
1875 pm2->Draw();
1876 return;
1877 }
1878 }
1879 }
1880 fTimer->Stop();
1881 fTimer->Print();
1882}
1883
1884////////////////////////////////////////////////////////////////////////////////
1885/// Check of validity of the normal for a given shape.
1886/// Sample points inside the shape. Generate directions randomly in cos(theta)
1887/// and propagate to boundary. Compute normal and safety at crossing point, plot
1888/// the point and generate a random direction so that (dir) dot (norm) <0.
1889
1891{
1892 Double_t dx = ((TGeoBBox*)shape)->GetDX();
1893 Double_t dy = ((TGeoBBox*)shape)->GetDY();
1894 Double_t dz = ((TGeoBBox*)shape)->GetDZ();
1895 Double_t dmax = 2.*TMath::Sqrt(dx*dx+dy*dy+dz*dz);
1896 // Number of tracks shot for every point inside the shape
1897 const Int_t kNtracks = 1000;
1898 Int_t n10 = nsamples/10;
1899 Int_t itot = 0, errcnt = 0, errsame=0;
1900 Int_t i;
1901 Double_t dist, olddist, safe, dot;
1902 Double_t theta, phi, ndotd;
1903 Double_t *spoint = new Double_t[3*nsamples];
1904 Double_t *sdir = new Double_t[3*nsamples];
1905 while (itot<nsamples) {
1906 Bool_t inside = kFALSE;
1907 while (!inside) {
1908 spoint[3*itot] = gRandom->Uniform(-dx,dx);
1909 spoint[3*itot+1] = gRandom->Uniform(-dy,dy);
1910 spoint[3*itot+2] = gRandom->Uniform(-dz,dz);
1911 inside = shape->Contains(&spoint[3*itot]);
1912 }
1913 phi = 2*TMath::Pi()*gRandom->Rndm();
1914 theta= TMath::ACos(1.-2.*gRandom->Rndm());
1915 sdir[3*itot] = TMath::Sin(theta)*TMath::Cos(phi);
1916 sdir[3*itot+1] = TMath::Sin(theta)*TMath::Sin(phi);
1917 sdir[3*itot+2] = TMath::Cos(theta);
1918 itot++;
1919 }
1920 Double_t point[3],newpoint[3], oldpoint[3];
1921 Double_t dir[3], olddir[3];
1922 Double_t norm[3], newnorm[3], oldnorm[3];
1923 TCanvas *errcanvas = 0;
1924 TPolyMarker3D *pm1 = 0;
1925 TPolyMarker3D *pm2 = 0;
1926 pm2 = new TPolyMarker3D();
1927// pm2->SetMarkerStyle(24);
1928 pm2->SetMarkerSize(0.2);
1929 pm2->SetMarkerColor(kBlue);
1930 if (!fTimer) fTimer = new TStopwatch();
1931 fTimer->Reset();
1932 fTimer->Start();
1933 for (itot = 0; itot<nsamples; itot++) {
1934 if (n10) {
1935 if ((itot%n10) == 0) printf("%i percent\n", Int_t(100*itot/nsamples));
1936 }
1937 oldnorm[0] = oldnorm[1] = oldnorm[2] = 0.;
1938 olddist = 0.;
1939 for (Int_t j=0; j<3; j++) {
1940 oldpoint[j] = point[j] = spoint[3*itot+j];
1941 olddir[j] = dir[j] = sdir[3*itot+j];
1942 }
1943 for (i=0; i<kNtracks; i++) {
1944 if (errcnt>0) break;
1945 dist = shape->DistFromInside(point,dir,3);
1946 for (Int_t j=0; j<3; j++) {
1947 newpoint[j] = point[j] + dist*dir[j];
1948 }
1949 shape->ComputeNormal(newpoint,dir,newnorm);
1950
1951 dot = olddir[0]*oldnorm[0]+olddir[1]*oldnorm[1]+ olddir[2]*oldnorm[2];
1952 if (!shape->Contains(point) && shape->Safety(point,kFALSE)>1.E-3) {
1953 errcnt++;
1954 printf("Error point outside (%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) =%g olddist=%g\n",
1955 point[0],point[1],point[2], dir[0], dir[1], dir[2], dist, olddist);
1956 printf(" old point: (%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f)\n",
1957 oldpoint[0],oldpoint[1],oldpoint[2], olddir[0], olddir[1], olddir[2]);
1958 if (!errcanvas) errcanvas = new TCanvas("shape_err03", Form("Shape %s (%s)",shape->GetName(),shape->ClassName()), 1000, 800);
1959 if (!pm1) {
1960 pm1 = new TPolyMarker3D();
1961 pm1->SetMarkerStyle(24);
1962 pm1->SetMarkerSize(0.4);
1963 pm1->SetMarkerColor(kRed);
1964 }
1965 pm1->SetNextPoint(point[0],point[1],point[2]);
1966 pm1->SetNextPoint(oldpoint[0],oldpoint[1],oldpoint[2]);
1967 break;
1968 }
1969 if ((dist<TGeoShape::Tolerance() && olddist*dot>1.E-3) || dist>dmax) {
1970 errsame++;
1971 if (errsame>1) {
1972 errcnt++;
1973 printf("Error DistFromInside(%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) =%g olddist=%g\n",
1974 point[0],point[1],point[2], dir[0], dir[1], dir[2], dist, olddist);
1975 printf(" new norm: (%g, %g, %g)\n", newnorm[0], newnorm[1], newnorm[2]);
1976 printf(" old point: (%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f)\n",
1977 oldpoint[0],oldpoint[1],oldpoint[2], olddir[0], olddir[1], olddir[2]);
1978 printf(" old norm: (%g, %g, %g)\n", oldnorm[0], oldnorm[1], oldnorm[2]);
1979 if (!errcanvas) errcanvas = new TCanvas("shape_err03", Form("Shape %s (%s)",shape->GetName(),shape->ClassName()), 1000, 800);
1980 if (!pm1) {
1981 pm1 = new TPolyMarker3D();
1982 pm1->SetMarkerStyle(24);
1983 pm1->SetMarkerSize(0.4);
1984 pm1->SetMarkerColor(kRed);
1985 }
1986 pm1->SetNextPoint(point[0],point[1],point[2]);
1987 pm1->SetNextPoint(oldpoint[0],oldpoint[1],oldpoint[2]);
1988 break;
1989 }
1990 } else errsame = 0;
1991
1992 olddist = dist;
1993 for (Int_t j=0; j<3; j++) {
1994 oldpoint[j] = point[j];
1995 point[j] += dist*dir[j];
1996 }
1997 safe = shape->Safety(point, kTRUE);
1998 if (safe>1.E-3) {
1999 errcnt++;
2000 printf("Error safety (%19.15f, %19.15f, %19.15f) safe=%g\n",
2001 point[0],point[1],point[2], safe);
2002 if (!errcanvas) errcanvas = new TCanvas("shape_err03", Form("Shape %s (%s)",shape->GetName(),shape->ClassName()), 1000, 800);
2003 if (!pm1) {
2004 pm1 = new TPolyMarker3D();
2005 pm1->SetMarkerStyle(24);
2006 pm1->SetMarkerSize(0.4);
2007 pm1->SetMarkerColor(kRed);
2008 }
2009 pm1->SetNextPoint(point[0],point[1],point[2]);
2010 break;
2011 }
2012 // Compute normal
2013 shape->ComputeNormal(point,dir,norm);
2014 memcpy(oldnorm, norm, 3*sizeof(Double_t));
2015 memcpy(olddir, dir, 3*sizeof(Double_t));
2016 while (1) {
2017 phi = 2*TMath::Pi()*gRandom->Rndm();
2018 theta= TMath::ACos(1.-2.*gRandom->Rndm());
2019 dir[0]=TMath::Sin(theta)*TMath::Cos(phi);
2020 dir[1]=TMath::Sin(theta)*TMath::Sin(phi);
2021 dir[2]=TMath::Cos(theta);
2022 ndotd = dir[0]*norm[0]+dir[1]*norm[1]+dir[2]*norm[2];
2023 if (ndotd<0) break; // backwards, still inside shape
2024 }
2025 if ((itot%10) == 0) pm2->SetNextPoint(point[0],point[1],point[2]);
2026 }
2027 }
2028 fTimer->Stop();
2029 fTimer->Print();
2030 if (errcanvas) {
2031 shape->Draw();
2032 pm1->Draw();
2033 }
2034
2035 new TCanvas("shape03", Form("Shape %s (%s)",shape->GetName(),shape->ClassName()), 1000, 800);
2036 pm2->Draw();
2037 delete [] spoint;
2038 delete [] sdir;
2039}
2040
2041////////////////////////////////////////////////////////////////////////////////
2042/// Generate a lego plot fot the top volume, according to option.
2043
2045 Int_t nphi, Double_t phimin, Double_t phimax,
2046 Double_t /*rmin*/, Double_t /*rmax*/, Option_t *option)
2047{
2048 TH2F *hist = new TH2F("lego", option, nphi, phimin, phimax, ntheta, themin, themax);
2049
2050 Double_t degrad = TMath::Pi()/180.;
2051 Double_t theta, phi, step, matprop, x;
2052 Double_t start[3];
2053 Double_t dir[3];
2054 TGeoNode *startnode, *endnode;
2055 Int_t i; // loop index for phi
2056 Int_t j; // loop index for theta
2057 Int_t ntot = ntheta * nphi;
2058 Int_t n10 = ntot/10;
2059 Int_t igen = 0, iloop=0;
2060 printf("=== Lego plot sph. => nrays=%i\n", ntot);
2061 for (i=1; i<=nphi; i++) {
2062 for (j=1; j<=ntheta; j++) {
2063 igen++;
2064 if (n10) {
2065 if ((igen%n10) == 0) printf("%i percent\n", Int_t(100*igen/ntot));
2066 }
2067 x = 0;
2068 theta = hist->GetYaxis()->GetBinCenter(j);
2069 phi = hist->GetXaxis()->GetBinCenter(i)+1E-3;
2070 start[0] = start[1] = start[2] = 1E-3;
2071 dir[0]=TMath::Sin(theta*degrad)*TMath::Cos(phi*degrad);
2072 dir[1]=TMath::Sin(theta*degrad)*TMath::Sin(phi*degrad);
2073 dir[2]=TMath::Cos(theta*degrad);
2074 fGeoManager->InitTrack(&start[0], &dir[0]);
2075 startnode = fGeoManager->GetCurrentNode();
2076 if (fGeoManager->IsOutside()) startnode=0;
2077 if (startnode) {
2078 matprop = startnode->GetVolume()->GetMaterial()->GetRadLen();
2079 } else {
2080 matprop = 0.;
2081 }
2083// fGeoManager->IsStepEntering();
2084 // find where we end-up
2085 endnode = fGeoManager->Step();
2086 step = fGeoManager->GetStep();
2087 while (step<1E10) {
2088 // now see if we can make an other step
2089 iloop=0;
2090 while (!fGeoManager->IsEntering()) {
2091 iloop++;
2092 fGeoManager->SetStep(1E-3);
2093 step += 1E-3;
2094 endnode = fGeoManager->Step();
2095 }
2096 if (iloop>1000) printf("%i steps\n", iloop);
2097 if (matprop>0) {
2098 x += step/matprop;
2099 }
2100 if (endnode==0 && step>1E10) break;
2101 // generate an extra step to cross boundary
2102 startnode = endnode;
2103 if (startnode) {
2104 matprop = startnode->GetVolume()->GetMaterial()->GetRadLen();
2105 } else {
2106 matprop = 0.;
2107 }
2108
2110 endnode = fGeoManager->Step();
2111 step = fGeoManager->GetStep();
2112 }
2113 hist->Fill(phi, theta, x);
2114 }
2115 }
2116 return hist;
2117}
2118
2119////////////////////////////////////////////////////////////////////////////////
2120/// Draw random points in the bounding box of a volume.
2121
2123{
2124 if (!vol) return;
2125 vol->VisibleDaughters(kTRUE);
2126 vol->Draw();
2127 TString opt = option;
2128 opt.ToLower();
2129 TObjArray *pm = new TObjArray(128);
2130 TPolyMarker3D *marker = 0;
2131 const TGeoShape *shape = vol->GetShape();
2132 TGeoBBox *box = (TGeoBBox *)shape;
2133 Double_t dx = box->GetDX();
2134 Double_t dy = box->GetDY();
2135 Double_t dz = box->GetDZ();
2136 Double_t ox = (box->GetOrigin())[0];
2137 Double_t oy = (box->GetOrigin())[1];
2138 Double_t oz = (box->GetOrigin())[2];
2139 Double_t *xyz = new Double_t[3];
2140 printf("Random box : %f, %f, %f\n", dx, dy, dz);
2141 TGeoNode *node = 0;
2142 printf("Start... %i points\n", npoints);
2143 Int_t i=0;
2144 Int_t igen=0;
2145 Int_t ic = 0;
2146 Int_t n10 = npoints/10;
2147 Double_t ratio=0;
2148 while (igen<npoints) {
2149 xyz[0] = ox-dx+2*dx*gRandom->Rndm();
2150 xyz[1] = oy-dy+2*dy*gRandom->Rndm();
2151 xyz[2] = oz-dz+2*dz*gRandom->Rndm();
2153 igen++;
2154 if (n10) {
2155 if ((igen%n10) == 0) printf("%i percent\n", Int_t(100*igen/npoints));
2156 }
2157 node = fGeoManager->FindNode();
2158 if (!node) continue;
2159 if (!node->IsOnScreen()) continue;
2160 // draw only points in overlapping/non-overlapping volumes
2161 if (opt.Contains("many") && !node->IsOverlapping()) continue;
2162 if (opt.Contains("only") && node->IsOverlapping()) continue;
2163 ic = node->GetColour();
2164 if ((ic<0) || (ic>=128)) ic = 1;
2165 marker = (TPolyMarker3D*)pm->At(ic);
2166 if (!marker) {
2167 marker = new TPolyMarker3D();
2168 marker->SetMarkerColor(ic);
2169// marker->SetMarkerStyle(8);
2170// marker->SetMarkerSize(0.4);
2171 pm->AddAt(marker, ic);
2172 }
2173 marker->SetNextPoint(xyz[0], xyz[1], xyz[2]);
2174 i++;
2175 }
2176 printf("Number of visible points : %i\n", i);
2177 if (igen) ratio = (Double_t)i/(Double_t)igen;
2178 printf("efficiency : %g\n", ratio);
2179 for (Int_t m=0; m<128; m++) {
2180 marker = (TPolyMarker3D*)pm->At(m);
2181 if (marker) marker->Draw("SAME");
2182 }
2184 printf("---Daughters of %s made invisible.\n", fGeoManager->GetTopVolume()->GetName());
2185 printf("---Make them visible with : gGeoManager->GetTopVolume()->VisibleDaughters();\n");
2186 delete pm;
2187 delete [] xyz;
2188}
2189
2190////////////////////////////////////////////////////////////////////////////////
2191/// Randomly shoot nrays from point (startx,starty,startz) and plot intersections
2192/// with surfaces for current top node.
2193
2194void TGeoChecker::RandomRays(Int_t nrays, Double_t startx, Double_t starty, Double_t startz, const char *target_vol, Bool_t check_norm)
2195{
2196 TObjArray *pm = new TObjArray(128);
2197 TString starget = target_vol;
2198 TPolyLine3D *line = 0;
2199 TPolyLine3D *normline = 0;
2201// vol->VisibleDaughters(kTRUE);
2202
2203 Bool_t random = kFALSE;
2204 if (nrays<=0) {
2205 nrays = 100000;
2206 random = kTRUE;
2207 }
2208 Double_t start[3];
2209 Double_t dir[3];
2210 const Double_t *point = fGeoManager->GetCurrentPoint();
2211 vol->Draw();
2212 printf("Start... %i rays\n", nrays);
2213 TGeoNode *startnode, *endnode;
2214 Bool_t vis1,vis2;
2215 Int_t i=0;
2216 Int_t ipoint, inull;
2217 Int_t itot=0;
2218 Int_t n10=nrays/10;
2219 Double_t theta,phi, step, normlen;
2220 Double_t ox = ((TGeoBBox*)vol->GetShape())->GetOrigin()[0];
2221 Double_t oy = ((TGeoBBox*)vol->GetShape())->GetOrigin()[1];
2222 Double_t oz = ((TGeoBBox*)vol->GetShape())->GetOrigin()[2];
2223 Double_t dx = ((TGeoBBox*)vol->GetShape())->GetDX();
2224 Double_t dy = ((TGeoBBox*)vol->GetShape())->GetDY();
2225 Double_t dz = ((TGeoBBox*)vol->GetShape())->GetDZ();
2226 normlen = TMath::Max(dx,dy);
2227 normlen = TMath::Max(normlen,dz);
2228 normlen *= 0.05;
2229 while (itot<nrays) {
2230 itot++;
2231 inull = 0;
2232 ipoint = 0;
2233 if (n10) {
2234 if ((itot%n10) == 0) printf("%i percent\n", Int_t(100*itot/nrays));
2235 }
2236 if (random) {
2237 start[0] = ox-dx+2*dx*gRandom->Rndm();
2238 start[1] = oy-dy+2*dy*gRandom->Rndm();
2239 start[2] = oz-dz+2*dz*gRandom->Rndm();
2240 } else {
2241 start[0] = startx;
2242 start[1] = starty;
2243 start[2] = startz;
2244 }
2245 phi = 2*TMath::Pi()*gRandom->Rndm();
2246 theta= TMath::ACos(1.-2.*gRandom->Rndm());
2247 dir[0]=TMath::Sin(theta)*TMath::Cos(phi);
2248 dir[1]=TMath::Sin(theta)*TMath::Sin(phi);
2249 dir[2]=TMath::Cos(theta);
2250 startnode = fGeoManager->InitTrack(start[0],start[1],start[2], dir[0],dir[1],dir[2]);
2251 line = 0;
2252 if (fGeoManager->IsOutside()) startnode=0;
2253 vis1 = kFALSE;
2254 if (target_vol) {
2255 if (startnode && starget==startnode->GetVolume()->GetName()) vis1 = kTRUE;
2256 } else {
2257 if (startnode && startnode->IsOnScreen()) vis1 = kTRUE;
2258 }
2259 if (vis1) {
2260 line = new TPolyLine3D(2);
2261 line->SetLineColor(startnode->GetVolume()->GetLineColor());
2262 line->SetPoint(ipoint++, start[0], start[1], start[2]);
2263 i++;
2264 pm->Add(line);
2265 }
2266 while ((endnode=fGeoManager->FindNextBoundaryAndStep())) {
2267 step = fGeoManager->GetStep();
2268 if (step<TGeoShape::Tolerance()) inull++;
2269 else inull = 0;
2270 if (inull>5) break;
2271 const Double_t *normal = 0;
2272 if (check_norm) {
2273 normal = fGeoManager->FindNormalFast();
2274 if (!normal) break;
2275 }
2276 vis2 = kFALSE;
2277 if (target_vol) {
2278 if (starget==endnode->GetVolume()->GetName()) vis2 = kTRUE;
2279 } else if (endnode->IsOnScreen()) vis2 = kTRUE;
2280 if (ipoint>0) {
2281 // old visible node had an entry point -> finish segment
2282 line->SetPoint(ipoint, point[0], point[1], point[2]);
2283 if (!vis2 && check_norm) {
2284 normline = new TPolyLine3D(2);
2285 normline->SetLineColor(kBlue);
2286 normline->SetLineWidth(1);
2287 normline->SetPoint(0, point[0], point[1], point[2]);
2288 normline->SetPoint(1, point[0]+normal[0]*normlen,
2289 point[1]+normal[1]*normlen,
2290 point[2]+normal[2]*normlen);
2291 pm->Add(normline);
2292 }
2293 ipoint = 0;
2294 line = 0;
2295 }
2296 if (vis2) {
2297 // create new segment
2298 line = new TPolyLine3D(2);
2299 line->SetLineColor(endnode->GetVolume()->GetLineColor());
2300 line->SetPoint(ipoint++, point[0], point[1], point[2]);
2301 i++;
2302 if (check_norm) {
2303 normline = new TPolyLine3D(2);
2304 normline->SetLineColor(kBlue);
2305 normline->SetLineWidth(2);
2306 normline->SetPoint(0, point[0], point[1], point[2]);
2307 normline->SetPoint(1, point[0]+normal[0]*normlen,
2308 point[1]+normal[1]*normlen,
2309 point[2]+normal[2]*normlen);
2310 }
2311 pm->Add(line);
2312 if (!random) pm->Add(normline);
2313 }
2314 }
2315 }
2316 // draw all segments
2317 for (Int_t m=0; m<pm->GetEntriesFast(); m++) {
2318 line = (TPolyLine3D*)pm->At(m);
2319 if (line) line->Draw("SAME");
2320 }
2321 printf("number of segments : %i\n", i);
2322// fGeoManager->GetTopVolume()->VisibleDaughters(kFALSE);
2323// printf("---Daughters of %s made invisible.\n", fGeoManager->GetTopVolume()->GetName());
2324// printf("---Make them visible with : gGeoManager->GetTopVolume()->VisibleDaughters();\n");
2325 delete pm;
2326}
2327
2328////////////////////////////////////////////////////////////////////////////////
2329/// shoot npoints randomly in a box of 1E-5 around current point.
2330/// return minimum distance to points outside
2331/// make sure that path to current node is updated
2332/// get the response of tgeo
2333
2335 const char* g3path)
2336{
2337 TGeoNode *node = fGeoManager->FindNode();
2338 TGeoNode *nodegeo = 0;
2339 TGeoNode *nodeg3 = 0;
2340 TGeoNode *solg3 = 0;
2341 if (!node) {dist=-1; return 0;}
2342 Bool_t hasg3 = kFALSE;
2343 if (strlen(g3path)) hasg3 = kTRUE;
2344 TString geopath = fGeoManager->GetPath();
2345 dist = 1E10;
2346 TString common = "";
2347 // cd to common path
2348 Double_t point[3];
2349 Double_t closest[3];
2350 TGeoNode *node1 = 0;
2351 TGeoNode *node_close = 0;
2352 dist = 1E10;
2353 Double_t dist1 = 0;
2354 // initialize size of random box to epsil
2355 Double_t eps[3];
2356 eps[0] = epsil; eps[1]=epsil; eps[2]=epsil;
2357 const Double_t *pointg = fGeoManager->GetCurrentPoint();
2358 if (hasg3) {
2359 TString spath = geopath;
2360 TString name = "";
2361 Int_t index=0;
2362 while (index>=0) {
2363 index = spath.Index("/", index+1);
2364 if (index>0) {
2365 name = spath(0, index);
2366 if (strstr(g3path, name.Data())) {
2367 common = name;
2368 continue;
2369 } else break;
2370 }
2371 }
2372 // if g3 response was given, cd to common path
2373 if (strlen(common.Data())) {
2374 while (strcmp(fGeoManager->GetPath(), common.Data()) && fGeoManager->GetLevel()) {
2375 nodegeo = fGeoManager->GetCurrentNode();
2376 fGeoManager->CdUp();
2377 }
2378 fGeoManager->cd(g3path);
2379 solg3 = fGeoManager->GetCurrentNode();
2380 while (strcmp(fGeoManager->GetPath(), common.Data()) && fGeoManager->GetLevel()) {
2381 nodeg3 = fGeoManager->GetCurrentNode();
2382 fGeoManager->CdUp();
2383 }
2384 if (!nodegeo) return 0;
2385 if (!nodeg3) return 0;
2386 fGeoManager->cd(common.Data());
2388 Double_t xyz[3], local[3];
2389 for (Int_t i=0; i<npoints; i++) {
2390 xyz[0] = point[0] - eps[0] + 2*eps[0]*gRandom->Rndm();
2391 xyz[1] = point[1] - eps[1] + 2*eps[1]*gRandom->Rndm();
2392 xyz[2] = point[2] - eps[2] + 2*eps[2]*gRandom->Rndm();
2393 nodeg3->MasterToLocal(&xyz[0], &local[0]);
2394 if (!nodeg3->GetVolume()->Contains(&local[0])) continue;
2395 dist1 = TMath::Sqrt((xyz[0]-point[0])*(xyz[0]-point[0])+
2396 (xyz[1]-point[1])*(xyz[1]-point[1])+(xyz[2]-point[2])*(xyz[2]-point[2]));
2397 if (dist1<dist) {
2398 // save node and closest point
2399 dist = dist1;
2400 node_close = solg3;
2401 // make the random box smaller
2402 eps[0] = TMath::Abs(point[0]-pointg[0]);
2403 eps[1] = TMath::Abs(point[1]-pointg[1]);
2404 eps[2] = TMath::Abs(point[2]-pointg[2]);
2405 }
2406 }
2407 }
2408 if (!node_close) dist = -1;
2409 return node_close;
2410 }
2411
2412 // save current point
2413 memcpy(&point[0], pointg, 3*sizeof(Double_t));
2414 for (Int_t i=0; i<npoints; i++) {
2415 // generate a random point in MARS
2416 fGeoManager->SetCurrentPoint(point[0] - eps[0] + 2*eps[0]*gRandom->Rndm(),
2417 point[1] - eps[1] + 2*eps[1]*gRandom->Rndm(),
2418 point[2] - eps[2] + 2*eps[2]*gRandom->Rndm());
2419 // check if new node is different from the old one
2420 if (node1!=node) {
2421 dist1 = TMath::Sqrt((point[0]-pointg[0])*(point[0]-pointg[0])+
2422 (point[1]-pointg[1])*(point[1]-pointg[1])+(point[2]-pointg[2])*(point[2]-pointg[2]));
2423 if (dist1<dist) {
2424 dist = dist1;
2425 node_close = node1;
2426 memcpy(&closest[0], pointg, 3*sizeof(Double_t));
2427 // make the random box smaller
2428 eps[0] = TMath::Abs(point[0]-pointg[0]);
2429 eps[1] = TMath::Abs(point[1]-pointg[1]);
2430 eps[2] = TMath::Abs(point[2]-pointg[2]);
2431 }
2432 }
2433 }
2434 // restore the original point and path
2435 fGeoManager->FindNode(point[0], point[1], point[2]); // really needed ?
2436 if (!node_close) dist=-1;
2437 return node_close;
2438}
2439
2440////////////////////////////////////////////////////////////////////////////////
2441/// Shoot one ray from start point with direction (dirx,diry,dirz). Fills input array
2442/// with points just after boundary crossings.
2443
2444Double_t *TGeoChecker::ShootRay(Double_t *start, Double_t dirx, Double_t diry, Double_t dirz, Double_t *array, Int_t &nelem, Int_t &dim, Double_t *endpoint) const
2445{
2446// Int_t array_dimension = 3*dim;
2447 nelem = 0;
2448 Int_t istep = 0;
2449 if (!dim) {
2450 printf("empty input array\n");
2451 return array;
2452 }
2453// fGeoManager->CdTop();
2454 const Double_t *point = fGeoManager->GetCurrentPoint();
2455 TGeoNode *endnode;
2456 Bool_t is_entering;
2457 Double_t step, forward;
2458 Double_t dir[3];
2459 dir[0] = dirx;
2460 dir[1] = diry;
2461 dir[2] = dirz;
2462 fGeoManager->InitTrack(start, &dir[0]);
2464// printf("Start : (%f,%f,%f)\n", point[0], point[1], point[2]);
2466 step = fGeoManager->GetStep();
2467// printf("---next : at step=%f\n", step);
2468 if (step>1E10) return array;
2469 endnode = fGeoManager->Step();
2470 is_entering = fGeoManager->IsEntering();
2471 while (step<1E10) {
2472 if (endpoint) {
2473 forward = dirx*(endpoint[0]-point[0])+diry*(endpoint[1]-point[1])+dirz*(endpoint[2]-point[2]);
2474 if (forward<1E-3) {
2475// printf("exit : Passed start point. nelem=%i\n", nelem);
2476 return array;
2477 }
2478 }
2479 if (is_entering) {
2480 if (nelem>=dim) {
2481 Double_t *temparray = new Double_t[3*(dim+20)];
2482 memcpy(temparray, array, 3*dim*sizeof(Double_t));
2483 delete [] array;
2484 array = temparray;
2485 dim += 20;
2486 }
2487 memcpy(&array[3*nelem], point, 3*sizeof(Double_t));
2488// printf("%i (%f, %f, %f) step=%f\n", nelem, point[0], point[1], point[2], step);
2489 nelem++;
2490 } else {
2491 if (endnode==0 && step>1E10) {
2492// printf("exit : NULL endnode. nelem=%i\n", nelem);
2493 return array;
2494 }
2495 if (!fGeoManager->IsEntering()) {
2496// if (startnode) printf("stepping %f from (%f, %f, %f) inside %s...\n", step,point[0], point[1], point[2], startnode->GetName());
2497// else printf("stepping %f from (%f, %f, %f) OUTSIDE...\n", step,point[0], point[1], point[2]);
2498 istep = 0;
2499 }
2500 while (!fGeoManager->IsEntering()) {
2501 istep++;
2502 if (istep>1E3) {
2503// Error("ShootRay", "more than 1000 steps. Step was %f", step);
2504 nelem = 0;
2505 return array;
2506 }
2507 fGeoManager->SetStep(1E-5);
2508 endnode = fGeoManager->Step();
2509 }
2510 if (istep>0) printf("%i steps\n", istep);
2511 if (nelem>=dim) {
2512 Double_t *temparray = new Double_t[3*(dim+20)];
2513 memcpy(temparray, array, 3*dim*sizeof(Double_t));
2514 delete [] array;
2515 array = temparray;
2516 dim += 20;
2517 }
2518 memcpy(&array[3*nelem], point, 3*sizeof(Double_t));
2519// if (endnode) printf("%i (%f, %f, %f) step=%f\n", nelem, point[0], point[1], point[2], step);
2520 nelem++;
2521 }
2523 step = fGeoManager->GetStep();
2524// printf("---next at step=%f\n", step);
2525 endnode = fGeoManager->Step();
2526 is_entering = fGeoManager->IsEntering();
2527 }
2528 return array;
2529// printf("exit : INFINITE step. nelem=%i\n", nelem);
2530}
2531
2532////////////////////////////////////////////////////////////////////////////////
2533/// Check time of finding "Where am I" for n points.
2534
2536{
2537 Bool_t recheck = !strcmp(option, "RECHECK");
2538 if (recheck) printf("RECHECK\n");
2539 const TGeoShape *shape = fGeoManager->GetTopVolume()->GetShape();
2540 Double_t dx = ((TGeoBBox*)shape)->GetDX();
2541 Double_t dy = ((TGeoBBox*)shape)->GetDY();
2542 Double_t dz = ((TGeoBBox*)shape)->GetDZ();
2543 Double_t ox = (((TGeoBBox*)shape)->GetOrigin())[0];
2544 Double_t oy = (((TGeoBBox*)shape)->GetOrigin())[1];
2545 Double_t oz = (((TGeoBBox*)shape)->GetOrigin())[2];
2546 Double_t *xyz = new Double_t[3*npoints];
2547 TStopwatch *timer = new TStopwatch();
2548 printf("Random box : %f, %f, %f\n", dx, dy, dz);
2549 timer->Start(kFALSE);
2550 Int_t i;
2551 for (i=0; i<npoints; i++) {
2552 xyz[3*i] = ox-dx+2*dx*gRandom->Rndm();
2553 xyz[3*i+1] = oy-dy+2*dy*gRandom->Rndm();
2554 xyz[3*i+2] = oz-dz+2*dz*gRandom->Rndm();
2555 }
2556 timer->Stop();
2557 printf("Generation time :\n");
2558 timer->Print();
2559 timer->Reset();
2560 TGeoNode *node, *node1;
2561 printf("Start... %i points\n", npoints);
2562 timer->Start(kFALSE);
2563 for (i=0; i<npoints; i++) {
2564 fGeoManager->SetCurrentPoint(xyz+3*i);
2565 if (recheck) fGeoManager->CdTop();
2566 node = fGeoManager->FindNode();
2567 if (recheck) {
2568 node1 = fGeoManager->FindNode();
2569 if (node1 != node) {
2570 printf("Difference for x=%g y=%g z=%g\n", xyz[3*i], xyz[3*i+1], xyz[3*i+2]);
2571 printf(" from top : %s\n", node->GetName());
2572 printf(" redo : %s\n", fGeoManager->GetPath());
2573 }
2574 }
2575 }
2576 timer->Stop();
2577 timer->Print();
2578 delete [] xyz;
2579 delete timer;
2580}
2581
2582////////////////////////////////////////////////////////////////////////////////
2583/// Geometry overlap checker based on sampling.
2584
2585void TGeoChecker::TestOverlaps(const char* path)
2586{
2588 printf("Checking overlaps for path :\n");
2589 if (!fGeoManager->cd(path)) return;
2590 TGeoNode *checked = fGeoManager->GetCurrentNode();
2591 checked->InspectNode();
2592 // shoot 1E4 points in the shape of the current volume
2593 Int_t npoints = 1000000;
2594 Double_t big = 1E6;
2595 Double_t xmin = big;
2596 Double_t xmax = -big;
2597 Double_t ymin = big;
2598 Double_t ymax = -big;
2599 Double_t zmin = big;
2600 Double_t zmax = -big;
2601 TObjArray *pm = new TObjArray(128);
2602 TPolyMarker3D *marker = 0;
2603 TPolyMarker3D *markthis = new TPolyMarker3D();
2604 markthis->SetMarkerColor(5);
2605 TNtuple *ntpl = new TNtuple("ntpl","random points","x:y:z");
2607 Double_t *point = new Double_t[3];
2608 Double_t dx = ((TGeoBBox*)shape)->GetDX();
2609 Double_t dy = ((TGeoBBox*)shape)->GetDY();
2610 Double_t dz = ((TGeoBBox*)shape)->GetDZ();
2611 Double_t ox = (((TGeoBBox*)shape)->GetOrigin())[0];
2612 Double_t oy = (((TGeoBBox*)shape)->GetOrigin())[1];
2613 Double_t oz = (((TGeoBBox*)shape)->GetOrigin())[2];
2614 Double_t *xyz = new Double_t[3*npoints];
2615 Int_t i=0;
2616 printf("Generating %i points inside %s\n", npoints, fGeoManager->GetPath());
2617 while (i<npoints) {
2618 point[0] = ox-dx+2*dx*gRandom->Rndm();
2619 point[1] = oy-dy+2*dy*gRandom->Rndm();
2620 point[2] = oz-dz+2*dz*gRandom->Rndm();
2621 if (!shape->Contains(point)) continue;
2622 // convert each point to MARS
2623// printf("local %9.3f %9.3f %9.3f\n", point[0], point[1], point[2]);
2624 fGeoManager->LocalToMaster(point, &xyz[3*i]);
2625// printf("master %9.3f %9.3f %9.3f\n", xyz[3*i], xyz[3*i+1], xyz[3*i+2]);
2626 xmin = TMath::Min(xmin, xyz[3*i]);
2627 xmax = TMath::Max(xmax, xyz[3*i]);
2628 ymin = TMath::Min(ymin, xyz[3*i+1]);
2629 ymax = TMath::Max(ymax, xyz[3*i+1]);
2630 zmin = TMath::Min(zmin, xyz[3*i+2]);
2631 zmax = TMath::Max(zmax, xyz[3*i+2]);
2632 i++;
2633 }
2634 delete [] point;
2635 ntpl->Fill(xmin,ymin,zmin);
2636 ntpl->Fill(xmax,ymin,zmin);
2637 ntpl->Fill(xmin,ymax,zmin);
2638 ntpl->Fill(xmax,ymax,zmin);
2639 ntpl->Fill(xmin,ymin,zmax);
2640 ntpl->Fill(xmax,ymin,zmax);
2641 ntpl->Fill(xmin,ymax,zmax);
2642 ntpl->Fill(xmax,ymax,zmax);
2643 ntpl->Draw("z:y:x");
2644
2645 // shoot the poins in the geometry
2646 TGeoNode *node;
2647 TString cpath;
2648 Int_t ic=0;
2649 TObjArray *overlaps = new TObjArray();
2650 printf("using FindNode...\n");
2651 for (Int_t j=0; j<npoints; j++) {
2652 // always start from top level (testing only)
2653 fGeoManager->CdTop();
2654 fGeoManager->SetCurrentPoint(&xyz[3*j]);
2655 node = fGeoManager->FindNode();
2656 cpath = fGeoManager->GetPath();
2657 if (cpath.Contains(path)) {
2658 markthis->SetNextPoint(xyz[3*j], xyz[3*j+1], xyz[3*j+2]);
2659 continue;
2660 }
2661 // current point is found in an overlapping node
2662 if (!node) ic=128;
2663 else ic = node->GetVolume()->GetLineColor();
2664 if (ic >= 128) ic = 0;
2665 marker = (TPolyMarker3D*)pm->At(ic);
2666 if (!marker) {
2667 marker = new TPolyMarker3D();
2668 marker->SetMarkerColor(ic);
2669 pm->AddAt(marker, ic);
2670 }
2671 // draw the overlapping point
2672 marker->SetNextPoint(xyz[3*j], xyz[3*j+1], xyz[3*j+2]);
2673 if (node) {
2674 if (overlaps->IndexOf(node) < 0) overlaps->Add(node);
2675 }
2676 }
2677 // draw all overlapping points
2678// for (Int_t m=0; m<128; m++) {
2679// marker = (TPolyMarker3D*)pm->At(m);
2680// if (marker) marker->Draw("SAME");
2681// }
2682 markthis->Draw("SAME");
2683 if (gPad) gPad->Update();
2684 // display overlaps
2685 if (overlaps->GetEntriesFast()) {
2686 printf("list of overlapping nodes :\n");
2687 for (i=0; i<overlaps->GetEntriesFast(); i++) {
2688 node = (TGeoNode*)overlaps->At(i);
2689 if (node->IsOverlapping()) printf("%s MANY\n", node->GetName());
2690 else printf("%s ONLY\n", node->GetName());
2691 }
2692 } else printf("No overlaps\n");
2693 delete ntpl;
2694 delete pm;
2695 delete [] xyz;
2696 delete overlaps;
2697}
2698
2699////////////////////////////////////////////////////////////////////////////////
2700/// Estimate weight of top level volume with a precision SIGMA(W)/W
2701/// better than PRECISION. Option can be "v" - verbose (default).
2702
2704{
2705 TList *matlist = fGeoManager->GetListOfMaterials();
2706 Int_t nmat = matlist->GetSize();
2707 if (!nmat) return 0;
2708 Int_t *nin = new Int_t[nmat];
2709 memset(nin, 0, nmat*sizeof(Int_t));
2710 TString opt = option;
2711 opt.ToLower();
2712 Bool_t isverbose = opt.Contains("v");
2714 Double_t dx = box->GetDX();
2715 Double_t dy = box->GetDY();
2716 Double_t dz = box->GetDZ();
2717 Double_t ox = (box->GetOrigin())[0];
2718 Double_t oy = (box->GetOrigin())[1];
2719 Double_t oz = (box->GetOrigin())[2];
2720 Double_t x,y,z;
2721 TGeoNode *node;
2722 TGeoMaterial *mat;
2723 Double_t vbox = 0.000008*dx*dy*dz; // m3
2724 Bool_t end = kFALSE;
2725 Double_t weight=0, sigma, eps, dens;
2726 Double_t eps0=1.;
2727 Int_t indmat;
2728 Int_t igen=0;
2729 Int_t iin = 0;
2730 while (!end) {
2731 x = ox-dx+2*dx*gRandom->Rndm();
2732 y = oy-dy+2*dy*gRandom->Rndm();
2733 z = oz-dz+2*dz*gRandom->Rndm();
2734 node = fGeoManager->FindNode(x,y,z);
2735 igen++;
2736 if (!node) continue;
2737 mat = node->GetVolume()->GetMedium()->GetMaterial();
2738 indmat = mat->GetIndex();
2739 if (indmat<0) continue;
2740 nin[indmat]++;
2741 iin++;
2742 if ((iin%100000)==0 || igen>1E8) {
2743 weight = 0;
2744 sigma = 0;
2745 for (indmat=0; indmat<nmat; indmat++) {
2746 mat = (TGeoMaterial*)matlist->At(indmat);
2747 dens = mat->GetDensity(); // [g/cm3]
2748 if (dens<1E-2) dens=0;
2749 dens *= 1000.; // [kg/m3]
2750 weight += dens*Double_t(nin[indmat]);
2751 sigma += dens*dens*nin[indmat];
2752 }
2754 eps = sigma/weight;
2755 weight *= vbox/Double_t(igen);
2756 sigma *= vbox/Double_t(igen);
2757 if (eps<precision || igen>1E8) {
2758 if (isverbose) {
2759 printf("=== Weight of %s : %g +/- %g [kg]\n",
2760 fGeoManager->GetTopVolume()->GetName(), weight, sigma);
2761 }
2762 end = kTRUE;
2763 } else {
2764 if (isverbose && eps<0.5*eps0) {
2765 printf("%8dK: %14.7g kg %g %%\n",
2766 igen/1000, weight, eps*100);
2767 eps0 = eps;
2768 }
2769 }
2770 }
2771 }
2772 delete [] nin;
2773 return weight;
2774}
2775////////////////////////////////////////////////////////////////////////////////
2776/// count voxel timing
2777
2779{
2780 TStopwatch timer;
2781 Double_t time;
2782 TGeoShape *shape = vol->GetShape();
2783 TGeoNode *node;
2784 TGeoMatrix *matrix;
2785 Double_t *point;
2786 Double_t local[3];
2787 Int_t *checklist;
2788 Int_t ncheck;
2790 TGeoStateInfo &td = *nav->GetCache()->GetInfo();
2791 timer.Start();
2792 for (Int_t i=0; i<npoints; i++) {
2793 point = xyz + 3*i;
2794 if (!shape->Contains(point)) continue;
2795 checklist = voxels->GetCheckList(point, ncheck, td);
2796 if (!checklist) continue;
2797 if (!ncheck) continue;
2798 for (Int_t id=0; id<ncheck; id++) {
2799 node = vol->GetNode(checklist[id]);
2800 matrix = node->GetMatrix();
2801 matrix->MasterToLocal(point, &local[0]);
2802 if (node->GetVolume()->GetShape()->Contains(&local[0])) break;
2803 }
2804 }
2805 nav->GetCache()->ReleaseInfo();
2806 time = timer.CpuTime();
2807 return time;
2808}
2809
2810////////////////////////////////////////////////////////////////////////////////
2811/// Returns optimal voxelization type for volume vol.
2812/// - kFALSE - cartesian
2813/// - kTRUE - cylindrical
2814
2816{
2817 return kFALSE;
2818}
unsigned int fFlags
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
#define e(i)
Definition RSha256.hxx:103
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
int Int_t
Definition RtypesCore.h:45
char Char_t
Definition RtypesCore.h:37
long Long_t
Definition RtypesCore.h:54
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
long long Long64_t
Definition RtypesCore.h:80
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
@ kRed
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
@ kYellow
Definition Rtypes.h:66
@ kFullCircle
Definition TAttMarker.h:55
Option_t Option_t option
Option_t Option_t cindex
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 r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize id
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
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 nchar
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
Option_t Option_t TPoint TPoint percent
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoIdentity * gGeoIdentity
Definition TGeoMatrix.h:478
float xmin
int nentries
float ymin
float xmax
float ymax
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2467
R__EXTERN TStyle * gStyle
Definition TStyle.h:414
#define gPad
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:33
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:38
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:40
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:45
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition TAxis.cxx:478
Generic 3D primitive description class.
Definition TBuffer3D.h:18
UInt_t NbPols() const
Definition TBuffer3D.h:82
UInt_t NbPnts() const
Definition TBuffer3D.h:80
UInt_t NbSegs() const
Definition TBuffer3D.h:81
TObject * fID
Definition TBuffer3D.h:87
Bool_t SetRawSizes(UInt_t reqPnts, UInt_t reqPntsCapacity, UInt_t reqSegs, UInt_t reqSegsCapacity, UInt_t reqPols, UInt_t reqPolsCapacity)
Set kRaw tessellation section of buffer with supplied sizes.
Double_t * fPnts
Definition TBuffer3D.h:112
The Canvas class.
Definition TCanvas.h:23
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition TFile.h:51
Box class.
Definition TGeoBBox.h:18
Geometry checking package.
Definition TGeoChecker.h:38
Int_t PropagateInGeom(Double_t *, Double_t *)
Propagate from START along DIR from boundary to boundary until exiting geometry.
Int_t fNchecks
Selected node for overlap checking.
Definition TGeoChecker.h:51
void CheckShape(TGeoShape *shape, Int_t testNo, Int_t nsamples, Option_t *option)
Test for shape navigation methods.
Bool_t * fFlags
Array of timing per volume.
Definition TGeoChecker.h:48
TStopwatch * fTimer
Array of flags per volume.
Definition TGeoChecker.h:49
TGeoVolume * fVsafe
Definition TGeoChecker.h:42
void CheckOverlapsBySampling(TGeoVolume *vol, Double_t ovlp=0.1, Int_t npoints=1000000) const
Check illegal overlaps for volume VOL within a limit OVLP by sampling npoints inside the volume shape...
TGeoChecker()
Default constructor.
virtual ~TGeoChecker()
Destructor.
void ShapeNormal(TGeoShape *shape, Int_t nsamples, Option_t *option)
Check of validity of the normal for a given shape.
void CheckOverlaps(const TGeoVolume *vol, Double_t ovlp=0.1, Option_t *option="")
Check illegal overlaps for volume VOL within a limit OVLP.
void TestOverlaps(const char *path)
Geometry overlap checker based on sampling.
TGeoOverlap * MakeCheckOverlap(const char *name, TGeoVolume *vol1, TGeoVolume *vol2, TGeoMatrix *mat1, TGeoMatrix *mat2, Bool_t isovlp, Double_t ovlp)
Check if the 2 non-assembly volume candidates overlap/extrude. Returns overlap object.
Double_t * fVal2
Array of number of crossings per volume.
Definition TGeoChecker.h:47
TGeoManager * fGeoManager
Definition TGeoChecker.h:41
void ShapeDistances(TGeoShape *shape, Int_t nsamples, Option_t *option)
Test TGeoShape::DistFromInside/Outside.
Int_t NChecksPerVolume(TGeoVolume *vol)
Compute number of overlaps combinations to check per volume.
void ShapeSafety(TGeoShape *shape, Int_t nsamples, Option_t *option)
Check of validity of safe distance for a given shape.
void CheckGeometryFull(Bool_t checkoverlaps=kTRUE, Bool_t checkcrossings=kTRUE, Int_t nrays=10000, const Double_t *vertex=nullptr)
Geometry checking.
void Score(TGeoVolume *, Int_t, Double_t)
Score a hit for VOL.
Double_t * ShootRay(Double_t *start, Double_t dirx, Double_t diry, Double_t dirz, Double_t *array, Int_t &nelem, Int_t &dim, Double_t *enpoint=nullptr) const
Shoot one ray from start point with direction (dirx,diry,dirz).
void CleanPoints(Double_t *points, Int_t &numPoints) const
Number of points on mesh to be checked.
TGeoNode * SamplePoints(Int_t npoints, Double_t &dist, Double_t epsil, const char *g3path)
shoot npoints randomly in a box of 1E-5 around current point.
TGeoNode * fSelectedNode
Timer.
Definition TGeoChecker.h:50
Double_t Weight(Double_t precision=0.01, Option_t *option="v")
Estimate weight of top level volume with a precision SIGMA(W)/W better than PRECISION.
void OpProgress(const char *opname, Long64_t current, Long64_t size, TStopwatch *watch=nullptr, Bool_t last=kFALSE, Bool_t refresh=kFALSE, const char *msg="")
Print current operation progress.
Int_t fNmeshPoints
Number of checks for current volume.
Definition TGeoChecker.h:52
void CheckGeometry(Int_t nrays, Double_t startx, Double_t starty, Double_t startz) const
Shoot nrays with random directions from starting point (startx, starty, startz) in the reference fram...
virtual void CheckBoundaryReference(Int_t icheck=-1)
Check the boundary errors reference file created by CheckBoundaryErrors method.
void CheckPoint(Double_t x=0, Double_t y=0, Double_t z=0, Option_t *option="")
Draw point (x,y,z) over the picture of the daughters of the volume containing this point.
TBuffer3D * fBuff2
Definition TGeoChecker.h:44
void SetNmeshPoints(Int_t npoints=1000)
Set number of points to be generated on the shape outline when checking for overlaps.
void PrintOverlaps() const
Print the current list of overlaps held by the manager class.
void Test(Int_t npoints, Option_t *option)
Check time of finding "Where am I" for n points.
Double_t TimingPerVolume(TGeoVolume *)
Compute timing per "FindNextBoundary" + "Safety" call.
Bool_t TestVoxels(TGeoVolume *vol, Int_t npoints=1000000)
Returns optimal voxelization type for volume vol.
Double_t * fVal1
Definition TGeoChecker.h:46
void RandomRays(Int_t nrays, Double_t startx, Double_t starty, Double_t startz, const char *target_vol=nullptr, Bool_t check_norm=kFALSE)
Randomly shoot nrays from point (startx,starty,startz) and plot intersections with surfaces for curre...
Bool_t fFullCheck
Definition TGeoChecker.h:45
TH2F * LegoPlot(Int_t ntheta=60, Double_t themin=0., Double_t themax=180., Int_t nphi=90, Double_t phimin=0., Double_t phimax=360., Double_t rmin=0., Double_t rmax=9999999, Option_t *option="")
Generate a lego plot fot the top volume, according to option.
void RandomPoints(TGeoVolume *vol, Int_t npoints, Option_t *option)
Draw random points in the bounding box of a volume.
virtual void CheckBoundaryErrors(Int_t ntracks=1000000, Double_t radius=-1.)
Check pushes and pulls needed to cross the next boundary with respect to the position given by FindNe...
Double_t CheckVoxels(TGeoVolume *vol, TGeoVoxelFinder *voxels, Double_t *xyz, Int_t npoints)
count voxel timing
TBuffer3D * fBuff1
Definition TGeoChecker.h:43
Matrix class used for computing global transformations Should NOT be used for node definition.
Definition TGeoMatrix.h:421
void Multiply(const TGeoMatrix *right)
multiply to the right with an other transformation if right is identity matrix, just return
A geometry iterator.
Definition TGeoNode.h:245
void Reset(TGeoVolume *top=nullptr)
Resets the iterator for volume TOP.
const TGeoMatrix * GetCurrentMatrix() const
Returns global matrix for current node.
void SetTopName(const char *name)
Set the top name for path.
Int_t GetLevel() const
Definition TGeoNode.h:276
void GetPath(TString &path) const
Returns the path for the current node.
TGeoNode * GetNode(Int_t level) const
Returns current node at a given level.
void Skip()
Stop iterating the current branch.
The manager class for any TGeo geometry.
Definition TGeoManager.h:45
TGeoNode * GetMother(Int_t up=1) const
Double_t * FindNormalFast()
Computes fast normal to next crossed boundary, assuming that the current point is close enough to the...
TObjArray * GetListOfUVolumes() const
TObjArray * GetListOfOverlaps()
void CdUp()
Go one level up in geometry.
virtual Bool_t cd(const char *path="")
Browse the tree of nodes starting from fTopNode according to pathname.
void LocalToMaster(const Double_t *local, Double_t *master) const
void RestoreMasterVolume()
Restore the master volume of the geometry.
TGeoNode * FindNextDaughterBoundary(Double_t *point, Double_t *dir, Int_t &idaughter, Bool_t compmatrix=kFALSE)
Computes as fStep the distance to next daughter of the current volume.
TGeoNavigator * GetCurrentNavigator() const
Returns current navigator for the calling thread.
TGeoVolume * GetMasterVolume() const
TGeoNode * GetCurrentNode() const
void SetCurrentDirection(Double_t *dir)
void SetVisLevel(Int_t level=3)
set default level down to which visualization is performed
TGeoNode * FindNextBoundary(Double_t stepmax=TGeoShape::Big(), const char *path="", Bool_t frombdr=kFALSE)
Find distance to next boundary and store it in fStep.
TGeoNode * Step(Bool_t is_geom=kTRUE, Bool_t cross=kTRUE)
Make a rectilinear step of length fStep from current point (fPoint) on current direction (fDirection)...
TGeoVolume * GetVolume(const char *name) const
Search for a named volume. All trailing blanks stripped.
TGeoNode * FindNextBoundaryAndStep(Double_t stepmax=TGeoShape::Big(), Bool_t compsafe=kFALSE)
Compute distance to next boundary within STEPMAX.
void SetCurrentPoint(Double_t *point)
Double_t * FindNormal(Bool_t forward=kTRUE)
Computes normal vector to the next surface that will be or was already crossed when propagating on a ...
TGeoNode * FindNode(Bool_t safe_start=kTRUE)
Returns deepest node containing current point.
const Double_t * GetCurrentPoint() const
TGeoVolume * MakeSphere(const char *name, TGeoMedium *medium, Double_t rmin, Double_t rmax, Double_t themin=0, Double_t themax=180, Double_t phimin=0, Double_t phimax=360)
Make in one step a volume pointing to a sphere shape with given medium.
Bool_t IsOutside() const
TGeoNode * InitTrack(const Double_t *point, const Double_t *dir)
Initialize current point and current direction vector (normalized) in MARS.
Int_t GetLevel() const
Double_t GetStep() const
TGeoHMatrix * GetCurrentMatrix() const
TGeoNode * GetTopNode() const
Bool_t IsSameLocation(Double_t x, Double_t y, Double_t z, Bool_t change=kFALSE)
Checks if point (x,y,z) is still in the current node.
const char * GetPath() const
Get path to the current node in the form /node0/node1/...
Int_t AddOverlap(const TNamed *ovlp)
Add an illegal overlap/extrusion to the list.
static void SetVerboseLevel(Int_t vl)
Return current verbosity level (static function).
void SetStep(Double_t step)
TGeoVolume * GetCurrentVolume() const
static Int_t GetVerboseLevel()
Set verbosity level (static function).
void CdTop()
Make top level node the current node.
Double_t Safety(Bool_t inside=kFALSE)
Compute safe distance from the current point.
void MasterToLocal(const Double_t *master, Double_t *local) const
void ResetState()
Reset current state flags.
void CdDown(Int_t index)
Make a daughter of current node current.
TList * GetListOfMaterials() const
TGeoVolume * GetTopVolume() const
void SetTopVisible(Bool_t vis=kTRUE)
make top volume visible on screen
void CheckOverlaps(Double_t ovlp=0.1, Option_t *option="")
Check all geometry for illegal overlaps within a limit OVLP.
Bool_t IsEntering() const
Base class describing materials.
Int_t GetIndex()
Retrieve material index in the list of materials.
virtual Double_t GetRadLen() const
virtual Double_t GetDensity() const
Geometrical transformation package.
Definition TGeoMatrix.h:41
virtual void LocalToMasterVect(const Double_t *local, Double_t *master) const
convert a vector by multiplying its column vector (x, y, z, 1) to matrix inverse
virtual void MasterToLocal(const Double_t *master, Double_t *local) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix
virtual void LocalToMaster(const Double_t *local, Double_t *master) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
TGeoMaterial * GetMaterial() const
Definition TGeoMedium.h:52
Class providing navigation API for TGeo geometries.
TGeoNodeCache * GetCache() const
TGeoStateInfo * GetInfo()
Get next state info pointer.
void ReleaseInfo()
Release last used state info pointer.
A node represent a volume positioned inside another.They store links to both volumes and to the TGeoM...
Definition TGeoNode.h:41
Bool_t IsOverlapping() const
Definition TGeoNode.h:105
Bool_t IsOnScreen() const
check if this node is drawn. Assumes that this node is current
Definition TGeoNode.cxx:273
TGeoVolume * GetVolume() const
Definition TGeoNode.h:97
void SetOverlaps(Int_t *ovlp, Int_t novlp)
set the list of overlaps for this node (ovlp must be created with operator new)
Definition TGeoNode.cxx:653
virtual TGeoMatrix * GetMatrix() const =0
Int_t GetColour() const
Definition TGeoNode.h:88
Int_t * GetOverlaps(Int_t &novlp) const
Definition TGeoNode.h:96
virtual void MasterToLocal(const Double_t *master, Double_t *local) const
Convert the point coordinates from mother reference to local reference system.
Definition TGeoNode.cxx:518
void InspectNode() const
Inspect this node.
Definition TGeoNode.cxx:282
Base class describing geometry overlaps.
Definition TGeoOverlap.h:41
TPolyMarker3D * GetPolyMarker() const
Definition TGeoOverlap.h:72
void SetNextPoint(Double_t x, Double_t y, Double_t z)
Set next overlapping point.
void SetOverlap(Double_t ovlp)
Definition TGeoOverlap.h:94
virtual void PrintInfo() const
Print some info.
Double_t GetOverlap() const
Definition TGeoOverlap.h:77
Base abstract class for all shapes.
Definition TGeoShape.h:26
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)=0
static Double_t Big()
Definition TGeoShape.h:88
virtual Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const =0
virtual void GetMeshNumbers(Int_t &, Int_t &, Int_t &) const
Definition TGeoShape.h:125
virtual Bool_t GetPointsOnSegments(Int_t npoints, Double_t *array) const =0
static void SetTransform(TGeoMatrix *matrix)
Set current transformation matrix that applies to shape.
virtual Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const =0
virtual const char * GetName() const
Get the shape name.
virtual Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const =0
virtual Double_t Capacity() const =0
virtual Bool_t Contains(const Double_t *point) const =0
static Double_t Tolerance()
Definition TGeoShape.h:91
virtual void SetPoints(Double_t *points) const =0
virtual void Draw(Option_t *option="")
Draw this shape.
Tessellated solid class.
Class describing translations.
Definition TGeoMatrix.h:122
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:49
TGeoMedium * GetMedium() const
Definition TGeoVolume.h:174
Bool_t Contains(const Double_t *point) const
Definition TGeoVolume.h:109
virtual TGeoNode * AddNode(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat=nullptr, Option_t *option="")
Add a TGeoNode to the list of nodes.
TGeoMaterial * GetMaterial() const
Definition TGeoVolume.h:173
Int_t GetNdaughters() const
Definition TGeoVolume.h:351
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
TObjArray * GetNodes()
Definition TGeoVolume.h:168
void VisibleDaughters(Bool_t vis=kTRUE)
set visibility for daughters
void FindOverlaps() const
loop all nodes marked as overlaps and find overlapping brothers
TGeoNode * GetNode(const char *name) const
get the pointer to a daughter node
virtual void SetVisibility(Bool_t vis=kTRUE)
set visibility of this volume
Int_t GetIndex(const TGeoNode *node) const
get index number for a given daughter
TGeoPatternFinder * GetFinder() const
Definition TGeoVolume.h:176
TGeoVoxelFinder * GetVoxels() const
Getter for optimization structure.
Int_t GetNumber() const
Definition TGeoVolume.h:183
TGeoShape * GetShape() const
Definition TGeoVolume.h:189
virtual void Draw(Option_t *option="")
draw top volume according to option
virtual Int_t GetCurrentNodeIndex() const
Definition TGeoVolume.h:166
virtual void DrawOnly(Option_t *option="")
draw only this volume
virtual void SetLineColor(Color_t lcolor)
Set the line color.
virtual Bool_t IsAssembly() const
Returns true if the volume is an assembly or a scaled assembly.
virtual Bool_t IsVisible() const
Definition TGeoVolume.h:154
void InspectShape() const
Definition TGeoVolume.h:194
Finder class handling voxels.
virtual Int_t * GetCheckList(const Double_t *point, Int_t &nelem, TGeoStateInfo &td)
get the list of daughter indices for which point is inside their bbox
virtual void Voxelize(Option_t *option="")
Voxelize attached volume according to option If the volume is an assembly, make sure the bbox is comp...
virtual void FindOverlaps(Int_t inode) const
create the list of nodes for which the bboxes overlap with inode's bbox
Bool_t NeedRebuild() const
1-D histogram with a double per channel (see TH1 documentation)}
Definition TH1.h:620
1-D histogram with a float per channel (see TH1 documentation)}
Definition TH1.h:577
virtual void LabelsOption(Option_t *option="h", Option_t *axis="X")
Sort bins with labels or set option(s) to draw axis with labels.
Definition TH1.cxx:5346
static void AddDirectory(Bool_t add=kTRUE)
Sets the flag controlling the automatic add of histograms in memory.
Definition TH1.cxx:1267
TAxis * GetXaxis()
Definition TH1.h:322
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3338
TAxis * GetYaxis()
Definition TH1.h:323
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3060
virtual UInt_t SetCanExtend(UInt_t extendBitMask)
Make the histogram axes extendable / not extendable according to the bit mask returns the previous bi...
Definition TH1.cxx:6631
@ kAllAxes
Definition TH1.h:75
virtual void LabelsDeflate(Option_t *axis="X")
Reduce the number of bins for the axis passed in the option to the number of bins having a label.
Definition TH1.cxx:5209
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition TH1.cxx:8856
2-D histogram with a float per channel (see TH1 documentation)}
Definition TH2.h:257
Int_t Fill(Double_t) override
Invalid Fill method.
Definition TH2.cxx:347
A doubly linked list.
Definition TList.h:38
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:357
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
A simple TTree restricted to a list of float variables only.
Definition TNtuple.h:28
Int_t Fill() override
Fill a Ntuple with current values in fArgs.
Definition TNtuple.cxx:169
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
Int_t IndexOf(const TObject *obj) const override
void AddAt(TObject *obj, Int_t idx) override
Add object at position ids.
Int_t GetEntries() const override
Return the number of objects in array (i.e.
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
TObject * RemoveAt(Int_t idx) override
Remove object at index idx.
void Add(TObject *obj) override
Definition TObjArray.h:68
Mother of all ROOT objects.
Definition TObject.h:41
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:207
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:956
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:970
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:274
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:944
A 3-dimensional polyline.
Definition TPolyLine3D.h:33
virtual void SetPoint(Int_t point, Double_t x, Double_t y, Double_t z)
Set point n to x, y, z.
A 3D polymarker.
virtual Int_t GetN() const
virtual Int_t SetNextPoint(Double_t x, Double_t y, Double_t z)
Set point following LastPoint to x, y, z.
void Draw(Option_t *option="") override
Draws 3-D polymarker with its current attributes.
Double_t Rndm() override
Machine independent random number generator.
Definition TRandom.cxx:552
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition TRandom.cxx:672
Stopwatch class.
Definition TStopwatch.h:28
Double_t RealTime()
Stop the stopwatch (if it is running) and return the realtime (in seconds) passed between the start a...
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Double_t CpuTime()
Stop the stopwatch (if it is running) and return the cputime (in seconds) passed between the start an...
void Continue()
Resume a stopped stopwatch.
void Stop()
Stop the stopwatch.
void Reset()
Definition TStopwatch.h:52
void Print(Option_t *option="") const override
Print the real and cpu time passed between the start and stop events.
Basic string class.
Definition TString.h:139
void ToLower()
Change string to lower-case.
Definition TString.cxx:1170
TSubString Strip(EStripType s=kTrailing, char c=' ') const
Return a substring of self stripped at beginning and/or end.
Definition TString.cxx:1151
const char * Data() const
Definition TString.h:380
@ kLeading
Definition TString.h:278
TString & Remove(Ssiz_t pos)
Definition TString.h:685
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2356
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:636
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition TString.h:651
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition TStyle.cxx:1589
A TTree represents a columnar dataset.
Definition TTree.h:79
virtual Int_t Fill()
Fill all branches.
Definition TTree.cxx:4594
virtual Int_t GetEntry(Long64_t entry, Int_t getall=0)
Read all branches of entry and return total number of bytes read.
Definition TTree.cxx:5629
void Draw(Option_t *opt) override
Default Draw method for all objects.
Definition TTree.h:428
virtual Int_t SetBranchAddress(const char *bname, void *add, TBranch **ptr=nullptr)
Change branch address, dealing with clone trees properly.
Definition TTree.cxx:8371
virtual Long64_t GetEntries() const
Definition TTree.h:460
TBranch * Branch(const char *name, T *obj, Int_t bufsize=32000, Int_t splitlevel=99)
Add a new branch, and infer the data type from the type of obj being passed.
Definition TTree.h:350
Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0) override
Write this object to the current directory.
Definition TTree.cxx:9708
TPaveText * pt
TLine * line
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition fillpatterns.C:1
const Double_t sigma
Double_t y[n]
Definition legend1.C:17
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
TH1F * h1
Definition legend1.C:5
return c2
Definition legend2.C:14
return c3
Definition legend3.C:15
Double_t ACos(Double_t)
Returns the principal value of the arc cosine of x, expressed in radians.
Definition TMath.h:630
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
Double_t Log(Double_t x)
Returns the natural logarithm of x.
Definition TMath.h:754
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:660
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:592
constexpr Double_t Pi()
Definition TMath.h:37
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:586
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
constexpr Double_t TwoPi()
Definition TMath.h:44
Statefull info for the current geometry level.
TMarker m
Definition textangle.C:8
#define dnext(otri1, otri2)
Definition triangle.c:987
#define lnext(otri1, otri2)
Definition triangle.c:943