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#### TGeoChecker::CheckOverlaps(Double_t ovlp, Option_t *option)
63Checks the geometry definition for illegal overlaps and illegal extrusions
64within a selected geometry hierarchy.
65
66An **illegal overlap** is reported when two placed volumes (typically sibling
67daughters of the same mother, or a daughter against the mother depending on
68the traversal) occupy a common region of space larger than the allowed tolerance.
69In practical terms, a set of test points generated on the surface of one
70candidate is found **inside** the other candidate volume by more than the configured
71tolerance.
72
73An **illegal extrusion** is reported when a placed daughter volume has surface
74points that lie outside its expected container region (e.g. outside its mother
75or outside the allowed placement region as defined by the geometry navigation
76context), again beyond the configured tolerance.
77
78The overlap checking is performed in three stages:
79
801. Candidate search and filtering
81The checker traverses the requested hierarchy and generates candidate pairs.
82Candidate pairs are aggressively reduced using fast bounding checks (including
83oriented bounding box filters) before any expensive navigation tests are performed.
84
852. Surface point generation and caching
86For each distinct TGeoShape appearing in the candidate list, the checker generates
87a set of surface sampling points and caches them. The cached sampling depends on
88the current meshing policy (see tuning below).
89
903. Exact overlap/extrusion tests (navigation-based)
91The final test uses TGeo navigation (Contains, Safety, frame transforms) to evaluate
92whether sampled surface points violate the overlap/extrusion constraints. This final
93stage can run in parallel when ROOT Implicit MT is enabled.
94
95* Usage:
96- Check the full geometry
97 - Call from the manager:
98
99~~~{.cpp}
100gGeoManager->CheckOverlaps(ovlp);
101~~~
102
103- Check only a sub-hierarchy
104 - Call from a volume:
105
106~~~{.cpp}
107gGeoManager->GetTopVolume()->CheckOverlaps(ovlp);
108~~~
109
110This checks overlaps/extrusions for the selected volume and its daughters
111(recursively), without requiring a full-geometry check.
112
113* Parameters:
114~~~{.cpp}
115ovlp
116~~~
117The allowed overlap tolerance. Violations larger than ovlp are reported
118as illegal overlaps/extrusions. The default value is 1.e-1, but a thorough overlap
119check may use 1e-6
120
121(The option string is intentionally not documented here; legacy sampling modes are
122deprecated in favor of dedicated sampling entry points.)
123
124* Tuning the surface sampling
125
126The quality and cost of the surface sampling used in overlap checking can be tuned
127via TGeoManager:
128
129~~~{.cpp}
130TGeoManager::SetNsegments(Int_t nseg) (default 20)
131~~~
132Controls the segmentation used by shapes whose surface discretization depends on
133angular subdivision (e.g. tubes, cones, polycones, polygonal shapes). Increasing
134nseg increases the geometric fidelity of the discretization and usually increases
135the time spent in the point-generation and checking stages.
136
137~~~{.cpp}
138TGeoManager::SetNmeshPoints(Int_t npoints) (default 1000)
139~~~
140Controls the number of additional surface points requested via GetPointsOnSegments.
141Increasing npoints improves sampling density (potentially catching smaller overlaps)
142at the cost of more computation in the point-generation and overlap checking stages.
143
144Both nseg and npoints affect the cached mesh points. Changing them changes the
145sampling policy and triggers regeneration of cached surface points at the next
146CheckOverlaps() call.
147
148* Parallel execution
149
150When ROOT Implicit MT is enabled via:
151~~~{.cpp}
152ROOT::EnableImplicitMT(N);
153~~~
154the checker runs the final navigation-based overlap/extrusion evaluation in parallel.
155Each worker thread initializes its own navigation state to avoid shared state in
156TGeo navigation.
157*/
158
159#include "TCanvas.h"
160#include "TStyle.h"
161#include "TFile.h"
162#include "TNtuple.h"
163#include "TH2.h"
164#include "TRandom3.h"
165#include "TPolyMarker3D.h"
166#include "TPolyLine3D.h"
167#include "TStopwatch.h"
168
169#include "TGeoVoxelFinder.h"
170#include "TGeoBBox.h"
171#include "TGeoPcon.h"
172#include "TGeoTessellated.h"
173#include "TGeoManager.h"
174#include "TGeoOverlap.h"
175#include "TGeoChecker.h"
176
177#include "TBuffer3D.h"
178#include "TBuffer3DTypes.h"
179#include "TMath.h"
180
181#include <cstdlib>
182
183// statics and globals
184
185////////////////////////////////////////////////////////////////////////////////
186/// Default constructor
187
190 fGeoManager(nullptr),
191 fVsafe(nullptr),
192 fFullCheck(kFALSE),
193 fVal1(nullptr),
194 fVal2(nullptr),
195 fFlags(nullptr),
196 fTimer(nullptr),
197 fSelectedNode(nullptr),
198 fNchecks(0)
199{
200}
201
202////////////////////////////////////////////////////////////////////////////////
203/// Constructor for a given geometry
204
207 fGeoManager(geom),
208 fVsafe(nullptr),
209 fFullCheck(kFALSE),
210 fVal1(nullptr),
211 fVal2(nullptr),
212 fFlags(nullptr),
213 fTimer(nullptr),
214 fSelectedNode(nullptr),
215 fNchecks(0)
216{
217}
218
219////////////////////////////////////////////////////////////////////////////////
220/// Destructor
221
223{
224 if (fTimer)
225 delete fTimer;
226}
227
228////////////////////////////////////////////////////////////////////////////////
229/// Print current operation progress.
230
232 Bool_t refresh, const char *msg)
233{
234 static Long64_t icount = 0;
235 static TString oname;
236 static TString nname;
237 static Long64_t ocurrent = 0;
238 static Long64_t osize = 0;
239 static Int_t oseconds = 0;
240 static TStopwatch *owatch = nullptr;
241 static Bool_t oneoftwo = kFALSE;
242 static Int_t nrefresh = 0;
243 const char symbol[4] = {'=', '\\', '|', '/'};
244 char progress[11] = " ";
245 Int_t ichar = icount % 4;
246 TString message(msg);
247 message += " ";
248
249 if (!refresh) {
250 nrefresh = 0;
251 if (!size)
252 return;
253 owatch = watch;
254 oname = opname;
255 ocurrent = TMath::Abs(current);
257 if (ocurrent > osize)
258 ocurrent = osize;
259 } else {
260 nrefresh++;
261 if (!osize)
262 return;
263 }
264 icount++;
265 Double_t time = 0.;
266 Int_t hours = 0;
267 Int_t minutes = 0;
268 Int_t seconds = 0;
269 if (owatch && !last) {
270 owatch->Stop();
271 time = owatch->RealTime();
272 hours = (Int_t)(time / 3600.);
273 time -= 3600 * hours;
274 minutes = (Int_t)(time / 60.);
275 time -= 60 * minutes;
276 seconds = (Int_t)time;
277 if (refresh) {
278 if (oseconds == seconds) {
279 owatch->Continue();
280 return;
281 }
283 }
284 oseconds = seconds;
285 }
286 if (refresh && oneoftwo) {
287 nname = oname;
288 if (fNchecks <= nrefresh)
289 fNchecks = nrefresh + 1;
290 Int_t pctdone = (Int_t)(100. * nrefresh / fNchecks);
291 oname = TString::Format(" == %3d%% ==", pctdone);
292 }
293 Double_t percent = 100.0 * ocurrent / osize;
294 Int_t nchar = Int_t(percent / 10);
295 if (nchar > 10)
296 nchar = 10;
297 Int_t i;
298 for (i = 0; i < nchar; i++)
299 progress[i] = '=';
300 progress[nchar] = symbol[ichar];
301 for (i = nchar + 1; i < 10; i++)
302 progress[i] = ' ';
303 progress[10] = '\0';
304 oname += " ";
305 oname.Remove(20);
306 if (size < 10000)
307 fprintf(stderr, "%s [%10s] %4lld ", oname.Data(), progress, ocurrent);
308 else if (size < 100000)
309 fprintf(stderr, "%s [%10s] %5lld ", oname.Data(), progress, ocurrent);
310 else
311 fprintf(stderr, "%s [%10s] %7lld ", oname.Data(), progress, ocurrent);
312 if (time > 0.)
313 fprintf(stderr, "[%6.2f %%] TIME %.2d:%.2d:%.2d %s\r", percent, hours, minutes, seconds, message.Data());
314 else
315 fprintf(stderr, "[%6.2f %%] %s\r", percent, message.Data());
316 if (refresh && oneoftwo)
317 oname = nname;
318 if (owatch)
319 owatch->Continue();
320 if (last) {
321 icount = 0;
322 owatch = nullptr;
323 ocurrent = 0;
324 osize = 0;
325 oseconds = 0;
327 nrefresh = 0;
328 fprintf(stderr, "\n");
329 }
330}
331
332////////////////////////////////////////////////////////////////////////////////
333/// Check pushes and pulls needed to cross the next boundary with respect to the
334/// position given by FindNextBoundary. If radius is not mentioned the full bounding
335/// box will be sampled.
336
338{
340 Info("CheckBoundaryErrors", "Top volume is %s", tvol->GetName());
341 const TGeoShape *shape = tvol->GetShape();
342 TGeoBBox *box = (TGeoBBox *)shape;
343 Double_t dl[3];
344 Double_t ori[3];
345 Double_t xyz[3];
346 Double_t nxyz[3];
347 Double_t dir[3];
349 Char_t path[1024];
350 Char_t cdir[10];
351
352 // Tree part
353 TFile *f = new TFile("geobugs.root", "recreate");
354 TTree *bug = new TTree("bug", "Geometrical problems");
355 bug->Branch("pos", xyz, "xyz[3]/D");
356 bug->Branch("dir", dir, "dir[3]/D");
357 bug->Branch("push", &relp, "push/D");
358 bug->Branch("path", &path, "path/C");
359 bug->Branch("cdir", &cdir, "cdir/C");
360
361 dl[0] = box->GetDX();
362 dl[1] = box->GetDY();
363 dl[2] = box->GetDZ();
364 ori[0] = (box->GetOrigin())[0];
365 ori[1] = (box->GetOrigin())[1];
366 ori[2] = (box->GetOrigin())[2];
367 if (radius > 0)
368 dl[0] = dl[1] = dl[2] = radius;
369
371 TH1F *hnew = new TH1F("hnew", "Precision pushing", 30, -20., 10.);
372 TH1F *hold = new TH1F("hold", "Precision pulling", 30, -20., 10.);
373 TH2F *hplotS = new TH2F("hplotS", "Problematic points", 100, -dl[0], dl[0], 100, -dl[1], dl[1]);
374 gStyle->SetOptStat(111111);
375
376 TGeoNode *node = nullptr;
377 Long_t igen = 0;
378 Long_t itry = 0;
379 Long_t n100 = ntracks / 100;
380 Double_t rad = TMath::Sqrt(dl[0] * dl[0] + dl[1] * dl[1]);
381 printf("Random box : %f, %f, %f, %f, %f, %f\n", ori[0], ori[1], ori[2], dl[0], dl[1], dl[2]);
382 printf("Start... %i points\n", ntracks);
383 if (!fTimer)
384 fTimer = new TStopwatch();
385 fTimer->Reset();
386 fTimer->Start();
387 while (igen < ntracks) {
389 Double_t r = rad * gRandom->Rndm();
390 xyz[0] = ori[0] + r * TMath::Cos(phi1);
391 xyz[1] = ori[1] + r * TMath::Sin(phi1);
392 Double_t z = (1. - 2. * gRandom->Rndm());
393 xyz[2] = ori[2] + dl[2] * z * TMath::Abs(z);
394 ++itry;
396 node = fGeoManager->FindNode();
397 if (!node || node == fGeoManager->GetTopNode())
398 continue;
399 ++igen;
400 if (n100 && !(igen % n100))
401 OpProgress("Sampling progress:", igen, ntracks, fTimer);
402 Double_t cost = 1. - 2. * gRandom->Rndm();
403 Double_t sint = TMath::Sqrt((1. + cost) * (1. - cost));
404 Double_t phi = TMath::TwoPi() * gRandom->Rndm();
405 dir[0] = sint * TMath::Cos(phi);
406 dir[1] = sint * TMath::Sin(phi);
407 dir[2] = cost;
410 Double_t step = fGeoManager->GetStep();
411
412 relp = 1.e-21;
413 for (Int_t i = 0; i < 30; ++i) {
414 relp *= 10.;
415 for (Int_t j = 0; j < 3; ++j)
416 nxyz[j] = xyz[j] + step * (1. + relp) * dir[j];
417 if (!fGeoManager->IsSameLocation(nxyz[0], nxyz[1], nxyz[2])) {
418 hnew->Fill(i - 20.);
419 if (i > 15) {
421 strncpy(path, fGeoManager->GetPath(), 1024);
422 path[1023] = '\0';
423 Double_t dotp = norm[0] * dir[0] + norm[1] * dir[1] + norm[2] * dir[2];
424 printf("Forward error i=%d p=%5.4f %5.4f %5.4f s=%5.4f dot=%5.4f path=%s\n", i, xyz[0], xyz[1], xyz[2],
425 step, dotp, path);
426 hplotS->Fill(xyz[0], xyz[1], (Double_t)i);
427 strncpy(cdir, "Forward", 10);
428 bug->Fill();
429 }
430 break;
431 }
432 }
433
434 relp = -1.e-21;
435 for (Int_t i = 0; i < 30; ++i) {
436 relp *= 10.;
437 for (Int_t j = 0; j < 3; ++j)
438 nxyz[j] = xyz[j] + step * (1. + relp) * dir[j];
439 if (fGeoManager->IsSameLocation(nxyz[0], nxyz[1], nxyz[2])) {
440 hold->Fill(i - 20.);
441 if (i > 15) {
443 strncpy(path, fGeoManager->GetPath(), 1024);
444 path[1023] = '\0';
445 Double_t dotp = norm[0] * dir[0] + norm[1] * dir[1] + norm[2] * dir[2];
446 printf("Backward error i=%d p=%5.4f %5.4f %5.4f s=%5.4f dot=%5.4f path=%s\n", i, xyz[0], xyz[1], xyz[2],
447 step, dotp, path);
448 strncpy(cdir, "Backward", 10);
449 bug->Fill();
450 }
451 break;
452 }
453 }
454 }
455 fTimer->Stop();
456
457 if (itry)
458 printf("CPU time/point = %5.2emus: Real time/point = %5.2emus\n", 1000000. * fTimer->CpuTime() / itry,
459 1000000. * fTimer->RealTime() / itry);
460 bug->Write();
461 delete bug;
462 bug = nullptr;
463 delete f;
464
466
467 if (itry)
468 printf("Effic = %3.1f%%\n", (100. * igen) / itry);
469 TCanvas *c1 = new TCanvas("c1", "Results", 600, 800);
470 c1->Divide(1, 2);
471 c1->cd(1);
472 gPad->SetLogy();
473 hold->Draw();
474 c1->cd(2);
475 gPad->SetLogy();
476 hnew->Draw();
477 /*TCanvas *c3 = */ new TCanvas("c3", "Plot", 600, 600);
478 hplotS->Draw("cont0");
479}
480
481////////////////////////////////////////////////////////////////////////////////
482/// Check the boundary errors reference file created by CheckBoundaryErrors method.
483/// The shape for which the crossing failed is drawn with the starting point in red
484/// and the extrapolated point to boundary (+/- failing push/pull) in yellow.
485
487{
488 Double_t xyz[3];
489 Double_t nxyz[3];
490 Double_t dir[3];
491 Double_t lnext[3];
492 Double_t push;
493 Char_t path[1024];
494 Char_t cdir[10];
495 // Tree part
496 TFile *f = new TFile("geobugs.root", "read");
497 TTree *bug = (TTree *)f->Get("bug");
498 bug->SetBranchAddress("pos", xyz);
499 bug->SetBranchAddress("dir", dir);
500 bug->SetBranchAddress("push", &push);
501 bug->SetBranchAddress("path", &path);
502 bug->SetBranchAddress("cdir", &cdir);
503
504 Int_t nentries = (Int_t)bug->GetEntries();
505 printf("nentries %d\n", nentries);
506 if (icheck < 0) {
507 for (Int_t i = 0; i < nentries; i++) {
508 bug->GetEntry(i);
509 printf("%-9s error push=%g p=%5.4f %5.4f %5.4f s=%5.4f dot=%5.4f path=%s\n", cdir, push, xyz[0], xyz[1],
510 xyz[2], 1., 1., path);
511 }
512 } else {
513 if (icheck >= nentries)
514 return;
517 bug->GetEntry(icheck);
518 printf("%-9s error push=%g p=%5.4f %5.4f %5.4f s=%5.4f dot=%5.4f path=%s\n", cdir, push, xyz[0], xyz[1], xyz[2],
519 1., 1., path);
524 Double_t step = fGeoManager->GetStep();
525 for (Int_t j = 0; j < 3; j++)
526 nxyz[j] = xyz[j] + step * (1. + 0.1 * push) * dir[j];
528 printf("step=%g in: %s\n", step, fGeoManager->GetPath());
529 printf(" -> next = %s push=%g change=%d\n", next->GetName(), push, (UInt_t)change);
530 next->GetVolume()->InspectShape();
531 next->GetVolume()->DrawOnly();
532 if (next != fGeoManager->GetCurrentNode()) {
534 if (index1 >= 0)
536 }
539 pm->SetNextPoint(lnext[0], lnext[1], lnext[2]);
540 pm->SetMarkerStyle(2);
541 pm->SetMarkerSize(0.2);
542 pm->SetMarkerColor(kRed);
543 pm->Draw("SAME");
545 for (Int_t j = 0; j < 3; j++)
546 nxyz[j] = xyz[j] + step * dir[j];
548 pm1->SetNextPoint(lnext[0], lnext[1], lnext[2]);
549 pm1->SetMarkerStyle(2);
550 pm1->SetMarkerSize(0.2);
551 pm1->SetMarkerColor(kYellow);
552 pm1->Draw("SAME");
554 }
555 delete bug;
556 delete f;
557}
558
559////////////////////////////////////////////////////////////////////////////////
560/// Geometry checking. Optional overlap checkings (by sampling and by mesh). Optional
561/// boundary crossing check + timing per volume.
562///
563/// STAGE 1: extensive overlap checking by sampling per volume. Stdout need to be
564/// checked by user to get report, then TGeoVolume::CheckOverlaps(0.01, "s") can
565/// be called for the suspicious volumes.
566///
567/// STAGE 2: normal overlap checking using the shapes mesh - fills the list of
568/// overlaps.
569///
570/// STAGE 3: shooting NRAYS rays from VERTEX and counting the total number of
571/// crossings per volume (rays propagated from boundary to boundary until
572/// geometry exit). Timing computed and results stored in a histo.
573///
574/// STAGE 4: shooting 1 mil. random rays inside EACH volume and calling
575/// FindNextBoundary() + Safety() for each call. The timing is normalized by the
576/// number of crossings computed at stage 2 and presented as percentage.
577/// One can get a picture on which are the most "burned" volumes during
578/// transportation from geometry point of view. Another plot of the timing per
579/// volume vs. number of daughters is produced.
580///
581/// All histos are saved in the file statistics.root
582
584{
586 if (!fTimer)
587 fTimer = new TStopwatch();
588 Int_t i;
590 fFlags = new Bool_t[nuid];
591 memset(fFlags, 0, nuid * sizeof(Bool_t));
592 TGeoVolume *vol;
593 TCanvas *c = new TCanvas("overlaps", "Overlaps by sampling", 800, 800);
594
595 // STAGE 1: Overlap checking by sampling per volume
596 if (checkoverlaps) {
597 printf("====================================================================\n");
598 printf("STAGE 1: Overlap checking by sampling within 10 microns\n");
599 printf("====================================================================\n");
600 fGeoManager->CheckOverlaps(0.001, "s");
601
602 // STAGE 2: Global overlap/extrusion checking
603 printf("====================================================================\n");
604 printf("STAGE 2: Global overlap/extrusion checking within 10 microns\n");
605 printf("====================================================================\n");
607 }
608
609 if (!checkcrossings) {
610 delete[] fFlags;
611 fFlags = nullptr;
612 delete c;
613 return;
614 }
615
616 fVal1 = new Double_t[nuid];
617 fVal2 = new Double_t[nuid];
618 memset(fVal1, 0, nuid * sizeof(Double_t));
619 memset(fVal2, 0, nuid * sizeof(Double_t));
620 // STAGE 3: How many crossings per volume in a realistic case ?
621 // Ignore volumes with no daughters
622
623 // Generate rays from vertex in phi=[0,2*pi] theta=[0,pi]
624 // Int_t ntracks = 1000000;
625 printf("====================================================================\n");
626 printf("STAGE 3: Propagating %i tracks starting from vertex\n and counting number of boundary crossings...\n",
627 ntracks);
628 printf("====================================================================\n");
629 Int_t nbound = 0;
630 Double_t theta, phi;
631 Double_t point[3], dir[3];
632 memset(point, 0, 3 * sizeof(Double_t));
633 if (vertex)
634 memcpy(point, vertex, 3 * sizeof(Double_t));
635
636 fTimer->Reset();
637 fTimer->Start();
638 for (i = 0; i < ntracks; i++) {
639 phi = 2. * TMath::Pi() * gRandom->Rndm();
640 theta = TMath::ACos(1. - 2. * gRandom->Rndm());
641 dir[0] = TMath::Sin(theta) * TMath::Cos(phi);
642 dir[1] = TMath::Sin(theta) * TMath::Sin(phi);
643 dir[2] = TMath::Cos(theta);
644 if ((i % 100) == 0)
645 OpProgress("Transporting tracks", i, ntracks, fTimer);
646 // if ((i%1000)==0) printf("... remaining tracks %i\n", ntracks-i);
647 nbound += PropagateInGeom(point, dir);
648 }
649 if (!nbound) {
650 printf("No boundary crossed\n");
651 return;
652 }
653 fTimer->Stop();
654 Double_t time1 = fTimer->CpuTime() * 1.E6;
655 Double_t time2 = time1 / ntracks;
657 OpProgress("Transporting tracks", ntracks, ntracks, fTimer, kTRUE);
658 printf("Time for crossing %i boundaries: %g [ms]\n", nbound, time1);
659 printf("Time per track for full geometry traversal: %g [ms], per crossing: %g [ms]\n", time2, time3);
660
661 // STAGE 4: How much time per volume:
662
663 printf("====================================================================\n");
664 printf("STAGE 4: How much navigation time per volume per next+safety call\n");
665 printf("====================================================================\n");
667 TGeoNode *current;
668 TString path;
669 vol = fGeoManager->GetTopVolume();
670 memset(fFlags, 0, nuid * sizeof(Bool_t));
672 timer.Start();
673 i = 0;
674 char volname[30];
675 strncpy(volname, vol->GetName(), 15);
676 volname[15] = '\0';
677 OpProgress(volname, i++, nuid, &timer);
678 Score(vol, 1, TimingPerVolume(vol));
679 while ((current = next())) {
680 vol = current->GetVolume();
681 Int_t uid = vol->GetNumber();
682 if (fFlags[uid])
683 continue;
684 fFlags[uid] = kTRUE;
685 next.GetPath(path);
686 fGeoManager->cd(path.Data());
687 strncpy(volname, vol->GetName(), 15);
688 volname[15] = '\0';
689 OpProgress(volname, i++, nuid, &timer);
690 Score(vol, 1, TimingPerVolume(vol));
691 }
692 OpProgress("STAGE 4 completed", i, nuid, &timer, kTRUE);
693 // Draw some histos
695 TCanvas *c1 = new TCanvas("c2", "ncrossings", 10, 10, 900, 500);
696 c1->SetGrid();
697 c1->SetTopMargin(0.15);
698 TFile *f = new TFile("statistics.root", "RECREATE");
699 TH1F *h = new TH1F("h", "number of boundary crossings per volume", 3, 0, 3);
700 h->SetStats(false);
701 h->SetFillColor(38);
702 h->SetCanExtend(TH1::kAllAxes);
703
704 memset(fFlags, 0, nuid * sizeof(Bool_t));
705 for (i = 0; i < nuid; i++) {
706 vol = fGeoManager->GetVolume(i);
707 if (!vol->GetNdaughters())
708 continue;
709 time_tot_pertrack += fVal1[i] * fVal2[i];
710 h->Fill(vol->GetName(), (Int_t)fVal1[i]);
711 }
712 time_tot_pertrack /= ntracks;
713 h->LabelsDeflate();
714 h->LabelsOption(">", "X");
715 h->Draw();
716
717 TCanvas *c2 = new TCanvas("c3", "time spent per volume in navigation", 10, 10, 900, 500);
718 c2->SetGrid();
719 c2->SetTopMargin(0.15);
720 TH2F *h2 = new TH2F("h2", "time per FNB call vs. ndaughters", 100, 0, 100, 100, 0, 15);
721 h2->SetStats(false);
722 h2->SetMarkerStyle(2);
723 TH1F *h1 = new TH1F("h1", "percent of time spent per volume", 3, 0, 3);
724 h1->SetStats(false);
725 h1->SetFillColor(38);
727 for (i = 0; i < nuid; i++) {
728 vol = fGeoManager->GetVolume(i);
729 if (!vol || !vol->GetNdaughters())
730 continue;
731 value = fVal1[i] * fVal2[i] / ntracks / time_tot_pertrack;
732 h1->Fill(vol->GetName(), value);
733 h2->Fill(vol->GetNdaughters(), fVal2[i]);
734 }
735 h1->LabelsDeflate();
736 h1->LabelsOption(">", "X");
737 h1->Draw();
738 TCanvas *c3 = new TCanvas("c4", "timing vs. ndaughters", 10, 10, 900, 500);
739 c3->SetGrid();
740 c3->SetTopMargin(0.15);
741 h2->Draw();
742 f->Write();
743 delete[] fFlags;
744 fFlags = nullptr;
745 delete[] fVal1;
746 fVal1 = nullptr;
747 delete[] fVal2;
748 fVal2 = nullptr;
749 delete fTimer;
750 fTimer = nullptr;
751 delete c;
752}
753
754////////////////////////////////////////////////////////////////////////////////
755/// Propagate from START along DIR from boundary to boundary until exiting
756/// geometry. Fill array of hits.
757
759{
760 fGeoManager->InitTrack(start, dir);
761 TGeoNode *current = nullptr;
762 Int_t nzero = 0;
763 Int_t nhits = 0;
764 while (!fGeoManager->IsOutside()) {
765 if (nzero > 3) {
766 // Problems in trying to cross a boundary
767 printf("Error in trying to cross boundary of %s\n", current->GetName());
768 return nhits;
769 }
771 if (!current || fGeoManager->IsOutside())
772 return nhits;
773 Double_t step = fGeoManager->GetStep();
774 if (step < 2. * TGeoShape::Tolerance()) {
775 nzero++;
776 continue;
777 } else
778 nzero = 0;
779 // Generate the hit
780 nhits++;
781 TGeoVolume *vol = current->GetVolume();
782 Score(vol, 0, 1.);
783 Int_t iup = 1;
785 while (mother && mother->GetVolume()->IsAssembly()) {
786 Score(mother->GetVolume(), 0, 1.);
788 }
789 }
790 return nhits;
791}
792
793////////////////////////////////////////////////////////////////////////////////
794/// Score a hit for VOL
795
797{
798 Int_t uid = vol->GetNumber();
799 switch (ifield) {
800 case 0: fVal1[uid] += value; break;
801 case 1: fVal2[uid] += value;
802 }
803}
804
805////////////////////////////////////////////////////////////////////////////////
806/// Compute timing per "FindNextBoundary" + "Safety" call. Volume must be
807/// in the current path.
808
810{
811 fTimer->Reset();
812 const TGeoShape *shape = vol->GetShape();
813 TGeoBBox *box = (TGeoBBox *)shape;
814 Double_t dx = box->GetDX();
815 Double_t dy = box->GetDY();
816 Double_t dz = box->GetDZ();
817 Double_t ox = (box->GetOrigin())[0];
818 Double_t oy = (box->GetOrigin())[1];
819 Double_t oz = (box->GetOrigin())[2];
820 Double_t point[3], dir[3], lpt[3], ldir[3];
821 Double_t pstep = 0.;
823 Double_t theta, phi;
825 fTimer->Start();
826 for (Int_t i = 0; i < 1000000; i++) {
827 lpt[0] = ox - dx + 2 * dx * gRandom->Rndm();
828 lpt[1] = oy - dy + 2 * dy * gRandom->Rndm();
829 lpt[2] = oz - dz + 2 * dz * gRandom->Rndm();
831 fGeoManager->SetCurrentPoint(point[0], point[1], point[2]);
832 phi = 2 * TMath::Pi() * gRandom->Rndm();
833 theta = TMath::ACos(1. - 2. * gRandom->Rndm());
834 ldir[0] = TMath::Sin(theta) * TMath::Cos(phi);
835 ldir[1] = TMath::Sin(theta) * TMath::Sin(phi);
836 ldir[2] = TMath::Cos(theta);
841 // dist = TGeoShape::Big();
842 if (!vol->IsAssembly()) {
843 Bool_t inside = vol->Contains(lpt);
844 if (!inside) {
845 // dist = vol->GetShape()->DistFromOutside(lpt,ldir,3,pstep);
846 // if (dist>=pstep) continue;
847 } else {
848 vol->GetShape()->DistFromInside(lpt, ldir, 3, pstep);
849 }
850
851 if (!vol->GetNdaughters())
852 vol->GetShape()->Safety(lpt, inside);
853 }
854 if (vol->GetNdaughters()) {
857 }
858 }
859 fTimer->Stop();
861 Int_t uid = vol->GetNumber();
862 Int_t ncrossings = (Int_t)fVal1[uid];
863 if (!vol->GetNdaughters())
864 printf("Time for volume %s (shape=%s): %g [ms] ndaughters=%d ncross=%d\n", vol->GetName(),
866 else
867 printf("Time for volume %s (assemb=%d): %g [ms] ndaughters=%d ncross=%d\n", vol->GetName(), vol->IsAssembly(),
869 return time_per_track;
870}
871
872////////////////////////////////////////////////////////////////////////////////
873/// Shoot nrays with random directions from starting point (startx, starty, startz)
874/// in the reference frame of this volume. Track each ray until exiting geometry, then
875/// shoot backwards from exiting point and compare boundary crossing points.
876
878{
879 Int_t i, j;
880 Double_t start[3], end[3];
881 Double_t dir[3];
882 Double_t dummy[3];
883 Double_t eps = 0.;
884 Double_t *array1 = new Double_t[3 * 1000];
885 Double_t *array2 = new Double_t[3 * 1000];
886 TObjArray *pma = new TObjArray();
888 pm = new TPolyMarker3D();
889 pm->SetMarkerColor(2); // error > eps
890 pm->SetMarkerStyle(8);
891 pm->SetMarkerSize(0.4);
892 pma->AddAt(pm, 0);
893 pm = new TPolyMarker3D();
894 pm->SetMarkerColor(4); // point not found
895 pm->SetMarkerStyle(8);
896 pm->SetMarkerSize(0.4);
897 pma->AddAt(pm, 1);
898 pm = new TPolyMarker3D();
899 pm->SetMarkerColor(6); // extra point back
900 pm->SetMarkerStyle(8);
901 pm->SetMarkerSize(0.4);
902 pma->AddAt(pm, 2);
904 Int_t dim1 = 1000, dim2 = 1000;
905 if ((startx == 0) && (starty == 0) && (startz == 0))
906 eps = 1E-3;
907 start[0] = startx + eps;
908 start[1] = starty + eps;
909 start[2] = startz + eps;
910 Int_t n10 = nrays / 10;
911 Double_t theta, phi;
912 Double_t dw, dwmin, dx, dy, dz;
914 for (i = 0; i < nrays; i++) {
915 if (n10) {
916 if ((i % n10) == 0)
917 printf("%i percent\n", Int_t(100 * i / nrays));
918 }
919 phi = 2 * TMath::Pi() * gRandom->Rndm();
920 theta = TMath::ACos(1. - 2. * gRandom->Rndm());
921 dir[0] = TMath::Sin(theta) * TMath::Cos(phi);
922 dir[1] = TMath::Sin(theta) * TMath::Sin(phi);
923 dir[2] = TMath::Cos(theta);
924 // shoot direct ray
925 nelem1 = nelem2 = 0;
926 // printf("DIRECT %i\n", i);
927 array1 = ShootRay(&start[0], dir[0], dir[1], dir[2], array1, nelem1, dim1);
928 if (!nelem1)
929 continue;
930 // for (j=0; j<nelem1; j++) printf("%i : %f %f %f\n", j, array1[3*j], array1[3*j+1], array1[3*j+2]);
931 memcpy(&end[0], &array1[3 * (nelem1 - 1)], 3 * sizeof(Double_t));
932 // shoot ray backwards
933 // printf("BACK %i\n", i);
934 array2 = ShootRay(&end[0], -dir[0], -dir[1], -dir[2], array2, nelem2, dim2, &start[0]);
935 if (!nelem2) {
936 printf("#### NOTHING BACK ###########################\n");
937 for (j = 0; j < nelem1; j++) {
938 pm = (TPolyMarker3D *)pma->At(0);
939 pm->SetNextPoint(array1[3 * j], array1[3 * j + 1], array1[3 * j + 2]);
940 }
941 continue;
942 }
943 // printf("BACKWARDS\n");
944 Int_t k = nelem2 >> 1;
945 for (j = 0; j < k; j++) {
946 memcpy(&dummy[0], &array2[3 * j], 3 * sizeof(Double_t));
947 memcpy(&array2[3 * j], &array2[3 * (nelem2 - 1 - j)], 3 * sizeof(Double_t));
948 memcpy(&array2[3 * (nelem2 - 1 - j)], &dummy[0], 3 * sizeof(Double_t));
949 }
950 // for (j=0; j<nelem2; j++) printf("%i : %f ,%f ,%f \n", j, array2[3*j], array2[3*j+1], array2[3*j+2]);
951 if (nelem1 != nelem2)
952 printf("### DIFFERENT SIZES : nelem1=%i nelem2=%i ##########\n", nelem1, nelem2);
953 ist1 = ist2 = 0;
954 // check first match
955
956 dx = array1[3 * ist1] - array2[3 * ist2];
957 dy = array1[3 * ist1 + 1] - array2[3 * ist2 + 1];
958 dz = array1[3 * ist1 + 2] - array2[3 * ist2 + 2];
959 dw = dx * dir[0] + dy * dir[1] + dz * dir[2];
962 // printf("%i : %s (%g, %g, %g)\n", ist1, fGeoManager->GetPath(),
963 // array1[3*ist1], array1[3*ist1+1], array1[3*ist1+2]);
964 if (TMath::Abs(dw) < 1E-4) {
965 // printf(" matching %i (%g, %g, %g)\n", ist2, array2[3*ist2], array2[3*ist2+1], array2[3*ist2+2]);
966 ist2++;
967 } else {
968 printf("### NOT MATCHING %i f:(%f, %f, %f) b:(%f %f %f) DCLOSE=%f\n", ist2, array1[3 * ist1],
969 array1[3 * ist1 + 1], array1[3 * ist1 + 2], array2[3 * ist2], array2[3 * ist2 + 1],
970 array2[3 * ist2 + 2], dw);
971 pm = (TPolyMarker3D *)pma->At(0);
972 pm->SetNextPoint(array2[3 * ist2], array2[3 * ist2 + 1], array2[3 * ist2 + 2]);
973 if (dw < 0) {
974 // first boundary missed on way back
975 } else {
976 // first boundary different on way back
977 ist2++;
978 }
979 }
980
981 while ((ist1 < nelem1 - 1) && (ist2 < nelem2)) {
984 // printf("%i : %s (%g, %g, %g)\n", ist1+1, fGeoManager->GetPath(),
985 // array1[3*ist1+3], array1[3*ist1+4], array1[3*ist1+5]);
986
987 dx = array1[3 * ist1 + 3] - array1[3 * ist1];
988 dy = array1[3 * ist1 + 4] - array1[3 * ist1 + 1];
989 dz = array1[3 * ist1 + 5] - array1[3 * ist1 + 2];
990 // distance to next point
991 dwmin = dx * dir[0] + dy * dir[1] + dz * dir[2];
992 while (ist2 < nelem2) {
993 ifound = 0;
994 dx = array2[3 * ist2] - array1[3 * ist1];
995 dy = array2[3 * ist2 + 1] - array1[3 * ist1 + 1];
996 dz = array2[3 * ist2 + 2] - array1[3 * ist1 + 2];
997 dw = dx * dir[0] + dy * dir[1] + dz * dir[2];
998 if (TMath::Abs(dw - dwmin) < 1E-4) {
999 ist1++;
1000 ist2++;
1001 break;
1002 }
1003 if (dw < dwmin) {
1004 // point found on way back. Check if close enough to ist1+1
1005 ifound++;
1006 dw = dwmin - dw;
1007 if (dw < 1E-4) {
1008 // point is matching ist1+1
1009 // printf(" matching %i (%g, %g, %g) DCLOSE=%g\n", ist2, array2[3*ist2],
1010 // array2[3*ist2+1], array2[3*ist2+2], dw);
1011 ist2++;
1012 ist1++;
1013 break;
1014 } else {
1015 // extra boundary found on way back
1018 pm = (TPolyMarker3D *)pma->At(2);
1019 pm->SetNextPoint(array2[3 * ist2], array2[3 * ist2 + 1], array2[3 * ist2 + 2]);
1020 printf("### EXTRA BOUNDARY %i : %s found at DCLOSE=%f\n", ist2, fGeoManager->GetPath(), dw);
1021 ist2++;
1022 continue;
1023 }
1024 } else {
1025 if (!ifound) {
1026 // point ist1+1 not found on way back
1029 pm = (TPolyMarker3D *)pma->At(1);
1030 pm->SetNextPoint(array2[3 * ist1 + 3], array2[3 * ist1 + 4], array2[3 * ist1 + 5]);
1031 printf("### BOUNDARY MISSED BACK #########################\n");
1032 ist1++;
1033 break;
1034 } else {
1035 ist1++;
1036 break;
1037 }
1038 }
1039 }
1040 }
1041 }
1042 pm = (TPolyMarker3D *)pma->At(0);
1043 pm->Draw("SAME");
1044 pm = (TPolyMarker3D *)pma->At(1);
1045 pm->Draw("SAME");
1046 pm = (TPolyMarker3D *)pma->At(2);
1047 pm->Draw("SAME");
1048 if (gPad) {
1049 gPad->Modified();
1050 gPad->Update();
1051 }
1052 delete[] array1;
1053 delete[] array2;
1054 delete pma; // markers are drawn on the pad
1055}
1056
1057////////////////////////////////////////////////////////////////////////////////
1058/// Clean-up the mesh of pcon/pgon from useless points
1059
1061{
1062 Int_t ipoint = 0;
1063 Int_t j, k = 0;
1064 Double_t rsq;
1065 for (Int_t i = 0; i < numPoints; i++) {
1066 j = 3 * i;
1067 rsq = points[j] * points[j] + points[j + 1] * points[j + 1];
1068 if (rsq < 1.e-10)
1069 continue;
1070 points[k] = points[j];
1071 points[k + 1] = points[j + 1];
1072 points[k + 2] = points[j + 2];
1073 ipoint++;
1074 k = 3 * ipoint;
1075 }
1076 numPoints = ipoint;
1077}
1078
1079////////////////////////////////////////////////////////////////////////////////
1080/// Compute overlap/extrusion for given candidate using mesh points of the shapes.
1081
1083{
1084 out = TGeoOverlapResult{}; // reset
1085
1086 TGeoVolume *vol1 = c.fVol1;
1087 TGeoVolume *vol2 = c.fVol2;
1088 if (!vol1 || !vol2)
1089 return kFALSE;
1090
1091 // Keep legacy early outs
1092 if (vol1->IsAssembly() || vol2->IsAssembly())
1093 return kFALSE;
1094
1095 if (dynamic_cast<TGeoTessellated *>(vol1->GetShape()) || dynamic_cast<TGeoTessellated *>(vol2->GetShape()))
1096 return kFALSE;
1097
1098 TGeoShape *shape1 = vol1->GetShape();
1099 TGeoShape *shape2 = vol2->GetShape();
1100
1101 const auto &m1 = GetMeshPoints(shape1);
1102 const auto &m2 = GetMeshPoints(shape2);
1103 if (m1.fNPoints <= 0 || m2.fNPoints <= 0)
1104 return kFALSE;
1105
1106 const Double_t *points1 = m1.fPoints.data();
1107 const Double_t *points2 = m2.fPoints.data();
1108 const Int_t numPoints1 = m1.fNPoints;
1109 const Int_t numPoints2 = m2.fNPoints;
1110
1111 const TGeoHMatrix mat1 = c.fMat1;
1112 const TGeoHMatrix mat2 = c.fMat2;
1113
1115 const Double_t ovlp = c.fOvlp;
1116
1117 auto add_point = [&](const Double_t p[3]) {
1118 if (out.fPoints.size() < 100)
1119 out.fPoints.push_back({{p[0], p[1], p[2]}});
1120 };
1121
1122 auto start_result_if_needed = [&]() {
1123 if (out.fVol1)
1124 return; // already started
1125 out.fName = c.fName;
1126 out.fVol1 = vol1;
1127 out.fVol2 = vol2;
1128 out.fMat1 = mat1;
1129 out.fMat2 = mat2;
1130 out.fIsOverlap = c.fIsOverlap;
1131 out.fMaxOverlap = 0.0;
1132 };
1133
1134 // Scratch
1135 Double_t local[3], local1[3], point[3];
1137
1138 // ---------------------------
1139 // Extrusion case (c.fIsOverlap == false):
1140 // "Test vol2 extrudes vol1" and also the legacy “mother points inside daughter” pass.
1141 // ---------------------------
1142 if (!c.fIsOverlap) {
1143 Bool_t found = kFALSE;
1144
1145 // loop all points of the daughter (vol2)
1146 for (Int_t ip = 0; ip < numPoints2; ++ip) {
1147 memcpy(local, &points2[3 * ip], 3 * sizeof(Double_t));
1149 continue;
1150
1151 mat2.LocalToMaster(local, point);
1152 mat1.MasterToLocal(point, local);
1153
1154 Bool_t extrude = !shape1->Contains(local);
1155 if (extrude) {
1156 safety = shape1->Safety(local, kFALSE);
1157 if (safety < ovlp)
1158 extrude = kFALSE;
1159 }
1160
1161 if (extrude) {
1163 found = kTRUE;
1164 if (safety > out.fMaxOverlap)
1165 out.fMaxOverlap = safety;
1166 add_point(point);
1167 }
1168 }
1169
1170 // loop all points of the mother (vol1)
1171 for (Int_t ip = 0; ip < numPoints1; ++ip) {
1172 memcpy(local, &points1[3 * ip], 3 * sizeof(Double_t));
1173 if (local[0] < 1e-10 && local[1] < 1e-10)
1174 continue;
1175
1176 mat1.LocalToMaster(local, point);
1177 mat2.MasterToLocal(point, local1);
1178
1179 Bool_t extrude = shape2->Contains(local1);
1180 if (extrude) {
1181 safety = shape1->Safety(local, kTRUE);
1182 if (safety > 1E-6) {
1183 extrude = kFALSE;
1184 } else {
1185 safety = shape2->Safety(local1, kTRUE);
1186 if (safety < ovlp)
1187 extrude = kFALSE;
1188 }
1189 }
1190
1191 if (extrude) {
1193 found = kTRUE;
1194 if (safety > out.fMaxOverlap)
1195 out.fMaxOverlap = safety;
1196 add_point(point);
1197 }
1198 }
1199
1200 return found;
1201 }
1202
1203 // ---------------------------
1204 // Overlap case
1205 // ---------------------------
1206 Bool_t found = kFALSE;
1207
1208 // loop all points of first candidate (vol1) in vol2
1209 for (Int_t ip = 0; ip < numPoints1; ++ip) {
1210 memcpy(local, &points1[3 * ip], 3 * sizeof(Double_t));
1211 if (local[0] < 1e-10 && local[1] < 1e-10)
1212 continue;
1213
1214 mat1.LocalToMaster(local, point);
1215 mat2.MasterToLocal(point, local);
1216
1217 Bool_t overlap = shape2->Contains(local);
1218 if (overlap) {
1219 safety = shape2->Safety(local, kTRUE);
1220 if (safety < ovlp)
1221 overlap = kFALSE;
1222 }
1223
1224 if (overlap) {
1226 found = kTRUE;
1227 if (safety > out.fMaxOverlap)
1228 out.fMaxOverlap = safety;
1229 add_point(point);
1230 }
1231 }
1232
1233 // loop all points of second candidate (vol2) in vol1
1234 for (Int_t ip = 0; ip < numPoints2; ++ip) {
1235 memcpy(local, &points2[3 * ip], 3 * sizeof(Double_t));
1236 if (local[0] < 1e-10 && local[1] < 1e-10)
1237 continue;
1238
1239 mat2.LocalToMaster(local, point);
1240 mat1.MasterToLocal(point, local);
1241
1242 Bool_t overlap = shape1->Contains(local);
1243 if (overlap) {
1244 safety = shape1->Safety(local, kTRUE);
1245 if (safety < ovlp)
1246 overlap = kFALSE;
1247 }
1248
1249 if (overlap) {
1251 found = kTRUE;
1252 if (safety > out.fMaxOverlap)
1253 out.fMaxOverlap = safety;
1254 add_point(point);
1255 }
1256 }
1257
1258 return found;
1259}
1260
1261////////////////////////////////////////////////////////////////////////////////
1262/// Check illegal overlaps for volume VOL within a limit OVLP by sampling npoints
1263/// inside the volume shape.
1264
1266{
1267 Int_t nd = vol->GetNdaughters();
1268 if (nd < 2)
1269 return;
1271 if (!voxels)
1272 return;
1273 TGeoBBox *box = (TGeoBBox *)vol->GetShape();
1274 TGeoShape *shape;
1275 TGeoNode *node;
1276 Double_t dx = box->GetDX();
1277 Double_t dy = box->GetDY();
1278 Double_t dz = box->GetDZ();
1279 Double_t pt[3];
1280 Double_t local[3];
1281 Int_t *check_list = nullptr;
1282 Int_t ncheck = 0;
1283 const Double_t *orig = box->GetOrigin();
1284 Int_t ipoint = 0;
1285 Int_t itry = 0;
1286 Int_t iovlp = 0;
1287 Int_t id = 0, id0 = 0, id1 = 0;
1288 Bool_t in, incrt;
1289 Double_t safe;
1290 TString name1 = "";
1291 TString name2 = "";
1292 TGeoOverlap **flags = nullptr;
1293 TGeoNode *node1, *node2;
1294 Int_t novlps = 0;
1296 // Int_t tid = TGeoManager::ThreadId();
1298 TGeoStateInfo &td = *nav->GetCache()->GetInfo();
1299 while (ipoint < npoints) {
1300 // Shoot randomly in the bounding box.
1301 pt[0] = orig[0] - dx + 2. * dx * gRandom->Rndm();
1302 pt[1] = orig[1] - dy + 2. * dy * gRandom->Rndm();
1303 pt[2] = orig[2] - dz + 2. * dz * gRandom->Rndm();
1304 if (!vol->Contains(pt)) {
1305 itry++;
1306 if (itry > 10000 && !ipoint) {
1307 Error("CheckOverlapsBySampling", "No point inside volume!!! - aborting");
1308 break;
1309 }
1310 continue;
1311 }
1312 // Check if the point is inside one or more daughters
1313 in = kFALSE;
1314 ipoint++;
1315 check_list = voxels->GetCheckList(pt, ncheck, td);
1316 if (!check_list || ncheck < 2)
1317 continue;
1318 for (id = 0; id < ncheck; id++) {
1319 id0 = check_list[id];
1320 node = vol->GetNode(id0);
1321 // Ignore MANY's
1322 if (node->IsOverlapping())
1323 continue;
1324 node->GetMatrix()->MasterToLocal(pt, local);
1325 shape = node->GetVolume()->GetShape();
1326 incrt = shape->Contains(local);
1327 if (!incrt)
1328 continue;
1329 if (!in) {
1330 in = kTRUE;
1331 id1 = id0;
1332 continue;
1333 }
1334 // The point is inside 2 or more daughters, check safety
1335 safe = shape->Safety(local, kTRUE);
1336 if (safe < ovlp)
1337 continue;
1338 // We really have found an overlap -> store the point in a container
1339 iovlp++;
1340 if (!novlps) {
1341 if (flags)
1342 delete[] flags; // should never happen
1343 flags = new TGeoOverlap *[nd * nd];
1344 memset(flags, 0, nd * nd * sizeof(TGeoOverlap *));
1345 }
1346 TGeoOverlap *nodeovlp = flags[nd * id1 + id0];
1347 if (!nodeovlp) {
1348 novlps++;
1349 node1 = vol->GetNode(id1);
1350 name1 = node1->GetName();
1351 mat1 = node1->GetMatrix();
1352 Int_t cindex = node1->GetVolume()->GetCurrentNodeIndex();
1353 while (cindex >= 0) {
1354 node1 = node1->GetVolume()->GetNode(cindex);
1355 name1 += TString::Format("/%s", node1->GetName());
1356 mat1.Multiply(node1->GetMatrix());
1357 cindex = node1->GetVolume()->GetCurrentNodeIndex();
1358 }
1359 node2 = vol->GetNode(id0);
1360 name2 = node2->GetName();
1361 mat2 = node2->GetMatrix();
1362 cindex = node2->GetVolume()->GetCurrentNodeIndex();
1363 while (cindex >= 0) {
1364 node2 = node2->GetVolume()->GetNode(cindex);
1365 name2 += TString::Format("/%s", node2->GetName());
1366 mat2.Multiply(node2->GetMatrix());
1367 cindex = node2->GetVolume()->GetCurrentNodeIndex();
1368 }
1369 nodeovlp = new TGeoOverlap(
1370 TString::Format("Volume %s: node %s overlapping %s", vol->GetName(), name1.Data(), name2.Data()),
1371 node1->GetVolume(), node2->GetVolume(), &mat1, &mat2, kTRUE, safe);
1372 flags[nd * id1 + id0] = nodeovlp;
1374 }
1375 // Max 100 points per marker
1376 if (nodeovlp->GetPolyMarker()->GetN() < 100)
1377 nodeovlp->SetNextPoint(pt[0], pt[1], pt[2]);
1378 if (nodeovlp->GetOverlap() < safe)
1379 nodeovlp->SetOverlap(safe);
1380 }
1381 }
1382 nav->GetCache()->ReleaseInfo();
1383 if (flags)
1384 delete[] flags;
1385 if (!novlps)
1386 return;
1387 Double_t capacity = vol->GetShape()->Capacity();
1388 capacity *= Double_t(iovlp) / Double_t(npoints);
1389 Double_t err = 1. / TMath::Sqrt(Double_t(iovlp));
1390 Info("CheckOverlapsBySampling", "#Found %d overlaps adding-up to %g +/- %g [cm3] for daughters of %s", novlps,
1391 capacity, err * capacity, vol->GetName());
1392}
1393
1394/// @brief Materialize a TGeoOverlapResult into a TGeoOverlap
1395/// @param r The result to materialize
1397{
1398 auto *ov = new TGeoOverlap(r.fName.Data(), r.fVol1, r.fVol2, &r.fMat1, &r.fMat2, r.fIsOverlap, r.fMaxOverlap);
1399
1400 for (const auto &p : r.fPoints)
1401 ov->SetNextPoint(p[0], p[1], p[2]);
1402
1404}
1405
1406////////////////////////////////////////////////////////////////////////////////
1407/// Create a policy stamp for overlap checking
1408
1416
1417////////////////////////////////////////////////////////////////////////////////
1418/// Compute number of overlaps combinations to check per volume
1419
1421{
1422 if (vol->GetFinder())
1423 return 0;
1424 UInt_t nd = vol->GetNdaughters();
1425 if (!nd)
1426 return 0;
1427 Bool_t is_assembly = vol->IsAssembly();
1428 TGeoIterator next1(vol);
1429 TGeoIterator next2(vol);
1430 Int_t nchecks = 0;
1431 TGeoNode *node;
1432 UInt_t id;
1433 if (!is_assembly) {
1434 while ((node = next1())) {
1435 if (node->IsOverlapping()) {
1436 next1.Skip();
1437 continue;
1438 }
1439 if (!node->GetVolume()->IsAssembly()) {
1440 nchecks++;
1441 next1.Skip();
1442 }
1443 }
1444 }
1445 // now check if the daughters overlap with each other
1446 if (nd < 2)
1447 return nchecks;
1448 TGeoVoxelFinder *vox = vol->GetVoxels();
1449 if (!vox)
1450 return nchecks;
1451
1453 Int_t novlp;
1454 Int_t *ovlps;
1455 Int_t ko;
1456 UInt_t io;
1457 for (id = 0; id < nd; id++) {
1458 node01 = vol->GetNode(id);
1459 if (node01->IsOverlapping())
1460 continue;
1461 vox->FindOverlaps(id);
1462 ovlps = node01->GetOverlaps(novlp);
1463 if (!ovlps)
1464 continue;
1465 for (ko = 0; ko < novlp; ko++) { // loop all possible overlapping candidates
1466 io = ovlps[ko]; // (node1, shaped, matrix1, points, fBuff1)
1467 if (io <= id)
1468 continue;
1469 node02 = vol->GetNode(io);
1470 if (node02->IsOverlapping())
1471 continue;
1472
1473 // We have to check node against node1, but they may be assemblies
1474 if (node01->GetVolume()->IsAssembly()) {
1475 next1.Reset(node01->GetVolume());
1476 while ((node = next1())) {
1477 if (!node->GetVolume()->IsAssembly()) {
1478 if (node02->GetVolume()->IsAssembly()) {
1479 next2.Reset(node02->GetVolume());
1480 while ((node1 = next2())) {
1481 if (!node1->GetVolume()->IsAssembly()) {
1482 nchecks++;
1483 next2.Skip();
1484 }
1485 }
1486 } else {
1487 nchecks++;
1488 }
1489 next1.Skip();
1490 }
1491 }
1492 } else {
1493 // node not assembly
1494 if (node02->GetVolume()->IsAssembly()) {
1495 next2.Reset(node02->GetVolume());
1496 while ((node1 = next2())) {
1497 if (!node1->GetVolume()->IsAssembly()) {
1498 nchecks++;
1499 next2.Skip();
1500 }
1501 }
1502 } else {
1503 // node1 also not assembly
1504 nchecks++;
1505 }
1506 }
1507 }
1508 node01->SetOverlaps(nullptr, 0);
1509 }
1510 return nchecks;
1511}
1512
1514{
1516 if (!shape)
1517 return kEmpty;
1518
1519 std::shared_lock lock(fMeshCacheMutex);
1520 auto it = fMeshPointsCache.find(shape);
1521 if (it != fMeshPointsCache.end())
1522 return it->second;
1523
1524 // Strict mode: caller forgot to include this shape in BuildMeshPointsCache
1525 ::Error("TGeoChecker::GetMeshPoints",
1526 "Mesh cache miss for shape %s. BuildMeshPointsCache did not cover all candidates.", shape->ClassName());
1527 return kEmpty;
1528}
1529
1530void TGeoChecker::BuildMeshPointsCache(const std::vector<TGeoOverlapCandidate> &candidates)
1531{
1532 std::unique_lock lock(fMeshCacheMutex);
1533
1534 const auto nowStamp = CurrentMeshPolicyStamp();
1536
1537 if (policyChanged) {
1538 fMeshPointsCache.clear();
1539 fMeshPointsCache.reserve(2 * candidates.size());
1541 }
1542
1543 auto ensureShape = [&](const TGeoShape *shape) {
1544 if (!shape)
1545 return;
1546 if (fMeshPointsCache.find(shape) != fMeshPointsCache.end())
1547 return;
1549
1550 // IMPORTANT: meshN is the *only* correct size for SetPoints() buffer.
1551 Int_t meshN = 0, meshS = 0, meshP = 0;
1552 const_cast<TGeoShape *>(shape)->GetMeshNumbers(meshN, meshS, meshP);
1553
1554 // 1) Base mesh points via SetPoints (exact sizing!)
1555 if (meshN > 0) {
1556 entry.fPoints.resize(3 * meshN);
1557 const_cast<TGeoShape *>(shape)->SetPoints(entry.fPoints.data());
1558 entry.fNPoints = meshN;
1559 } else {
1560 entry.fNPoints = 0;
1561 }
1562
1563 // 2) Extra points on segments: fill into exactly-sized array, then append
1565 if (nSegWanted > 0) {
1566 std::vector<Double_t> seg(3 * nSegWanted);
1567 Bool_t usedSeg = const_cast<TGeoShape *>(shape)->GetPointsOnSegments(nSegWanted, seg.data());
1568 if (usedSeg) {
1569 const auto oldSize = entry.fPoints.size();
1570 entry.fPoints.resize(oldSize + seg.size());
1571 std::memcpy(entry.fPoints.data() + oldSize, seg.data(), seg.size() * sizeof(Double_t));
1572 entry.fNPoints += nSegWanted;
1573 }
1574 }
1575
1576 // 3) Hardening: ensure consistency between count and storage.
1577 // (No dedup; just guard against inconsistencies.)
1578 if ((Int_t)entry.fPoints.size() != 3 * entry.fNPoints) {
1579 // This should never happen with the logic above, but keep a guard.
1580 const Int_t n3 = entry.fPoints.size() / 3;
1581 entry.fNPoints = n3;
1582 entry.fPoints.resize(3 * n3);
1583 }
1584
1585 fMeshPointsCache.emplace(shape, std::move(entry));
1586 };
1587
1588 // Ensure all shapes referenced by *this* candidate set exist in the cache.
1589 for (const auto &c : candidates) {
1590 if (c.fVol1)
1591 ensureShape(c.fVol1->GetShape());
1592 if (c.fVol2)
1593 ensureShape(c.fVol2->GetShape());
1594 }
1595}
1596
1599{
1600 // Cache hit: reuse last filled content for this shape
1601 if (lastShape == shape && cachedN > 0) {
1602 points = buff.fPnts;
1603 return cachedN;
1604 }
1605
1606 // Refill: compute mesh numbers only now
1607 Int_t meshN = buff.NbPnts();
1608 Int_t meshS = buff.NbSegs();
1609 Int_t meshP = buff.NbPols();
1610 shape->GetMeshNumbers(meshN, meshS, meshP);
1611
1613 buff.SetRawSizes(nraw, 3 * nraw, 0, 0, 0, 0);
1614 points = buff.fPnts;
1615
1617
1618 Int_t filledN = 0;
1619 if (useSeg) {
1621 } else {
1622 shape->SetPoints(points);
1623 filledN = meshN;
1624 }
1625
1626 lastShape = shape;
1627 cachedN = filledN;
1628 return filledN;
1629}
1630
1631////////////////////////////////////////////////////////////////////////////////
1632/// Stage1: Enumerate all overlap candidates for volume VOL within a limit OVLP.
1633
1635 std::vector<TGeoOverlapCandidate> &out)
1636{
1637 if (vol->GetFinder())
1638 return 0;
1639 UInt_t nd = vol->GetNdaughters();
1640 if (!nd)
1641 return 0;
1642
1643 const Bool_t is_assembly = vol->IsAssembly();
1644
1645 TString opt(option);
1646 opt.ToLower();
1647 fFullCheck = opt.Contains("f");
1648
1649 Int_t ncand = 0;
1650 Int_t nextr = 0;
1651 // Bool_t conv = vol->GetShape()->IsConvex();
1652
1653 // ---- EXTRUSIONS (only for daughters of a non-assembly volume)
1654 if (!is_assembly) {
1656 TGeoNode *node = nullptr;
1657 TGeoNode *nodecheck = nullptr;
1658 TString path;
1659 Int_t level;
1660
1661 while ((node = next1())) {
1662 if (node->IsOverlapping()) {
1663 next1.Skip();
1664 continue;
1665 }
1666 if (!node->GetVolume()->IsAssembly()) {
1667 if (fSelectedNode) {
1668 // only overlaps of the selected node OR real daughters if selected is assembly
1669 if ((fSelectedNode != node) && (!fSelectedNode->GetVolume()->IsAssembly())) {
1670 next1.Skip();
1671 continue;
1672 }
1673 if (node != fSelectedNode) {
1674 level = next1.GetLevel();
1675 while ((nodecheck = next1.GetNode(level--))) {
1676 if (nodecheck == fSelectedNode)
1677 break;
1678 }
1679 if (!nodecheck) {
1680 next1.Skip();
1681 continue;
1682 }
1683 }
1684 }
1685
1686 next1.GetPath(path);
1687 ncand++;
1688 nextr++;
1689 PushCandidate(out, TString::Format("%s extruded by: %s", vol->GetName(), path.Data()), (TGeoVolume *)vol,
1690 node->GetVolume(), gGeoIdentity, next1.GetCurrentMatrix(), kFALSE, ovlp);
1691
1692 next1.Skip();
1693 }
1694 }
1695 }
1696
1697 // if (nextr > 0)
1698 // printf("extrusion cand for %s : %d (convex = %d)\n", vol->GetName(), nextr, conv);
1699
1700 // ---- OVERLAPS between daughters
1701 if (nd < 2)
1702 return ncand;
1703
1704 TGeoVoxelFinder *vox = vol->GetVoxels();
1705 if (!vox) {
1706 Warning("EnumerateOverlapCandidates", "Volume %s with %i daughters but not voxelized", vol->GetName(), nd);
1707 return ncand;
1708 }
1709
1712
1713 TGeoNode *node1 = nullptr;
1714 TGeoNode *node01 = nullptr;
1715 TGeoNode *node02 = nullptr;
1716 TGeoNode *nodecheck = nullptr;
1717
1719 TGeoMatrix const *tr1, *tr2;
1720 TGeoBBox const *box1, *box2;
1721 TString path, path1;
1722 UInt_t id, io;
1723 Int_t level;
1724
1725 for (id = 0; id < nd; id++) {
1726 node01 = vol->GetNode(id);
1727 if (node01->IsOverlapping())
1728 continue;
1729
1730 next1.SetTopName(node01->GetName());
1731 path = node01->GetName();
1732 //=== ANY <-> ANY ===//
1733 for (io = id + 1; io < nd; io++) {
1734 // check node01 against node02 daughters of volume
1735 node02 = vol->GetNode(io);
1736 if (node02->IsOverlapping())
1737 continue;
1738 box1 = (const TGeoBBox *)node01->GetVolume()->GetShape();
1739 tr1 = node01->GetMatrix();
1740 box2 = (const TGeoBBox *)node02->GetVolume()->GetShape();
1741 tr2 = node02->GetMatrix();
1742 // OBB check
1744 continue;
1745
1746 next2.SetTopName(node02->GetName());
1747 path1 = node02->GetName();
1748
1749 // We have to check node01 against node02, but they may be assemblies
1750 if (node01->GetVolume()->IsAssembly()) {
1751 // left node assembly - make a visitor
1752 next1.Reset(node01->GetVolume());
1753 TGeoNode *node = nullptr; // will iterate components of node01
1754 //=== ASSEMBLY/ANY <-> ANY ===//
1755 while ((node = next1())) {
1756 hmat1 = node01->GetMatrix();
1757 hmat1 *= *next1.GetCurrentMatrix();
1758 box1 = (const TGeoBBox *)node->GetVolume()->GetShape();
1759 tr1 = &hmat1;
1760 // OBB check
1762 // No intersection, skip the full left branch
1763 next1.Skip();
1764 continue;
1765 }
1766
1767 if (!node->GetVolume()->IsAssembly()) {
1768 //=== ASSEMBLY/REAL <-> ANY ===//
1769 // Current daughter of node01 is real, get its path and transformation
1770 next1.GetPath(path);
1771
1772 if (node02->GetVolume()->IsAssembly()) {
1773 next2.Reset(node02->GetVolume());
1774 //=== ASSEMBLY/REAL <-> ASSEMBLY/ANY ===//
1775 while ((node1 = next2())) {
1776 hmat2 = node02->GetMatrix();
1777 hmat2 *= *next2.GetCurrentMatrix();
1778 box2 = (const TGeoBBox *)node1->GetVolume()->GetShape();
1779 // OBB check
1781 // No intersection, skip the full right branch
1782 next2.Skip();
1783 continue;
1784 }
1785 if (!node1->GetVolume()->IsAssembly()) {
1786 //=== ASSEMBLY/REAL <-> ASSEMBLY/REAL ===//
1787 // Selected node skip logic
1788 if (fSelectedNode) {
1789 if ((fSelectedNode != node) && (fSelectedNode != node1) &&
1791 next2.Skip();
1792 continue;
1793 }
1794 if ((fSelectedNode != node1) && (fSelectedNode != node)) {
1795 level = next2.GetLevel();
1796 while ((nodecheck = next2.GetNode(level--))) {
1797 if (nodecheck == fSelectedNode)
1798 break;
1799 }
1800 if (node02 == fSelectedNode)
1801 nodecheck = node02;
1802 if (!nodecheck) {
1803 level = next1.GetLevel();
1804 while ((nodecheck = next1.GetNode(level--))) {
1805 if (nodecheck == fSelectedNode)
1806 break;
1807 }
1808 }
1809 if (node01 == fSelectedNode)
1810 nodecheck = node01;
1811 if (!nodecheck) {
1812 next2.Skip();
1813 continue;
1814 }
1815 }
1816 }
1817
1818 next2.GetPath(path1);
1819
1820 ncand++;
1821 PushCandidate(out,
1822 TString::Format("%s/%s overlapping %s/%s", vol->GetName(), path.Data(),
1823 vol->GetName(), path1.Data()),
1824 node->GetVolume(), node1->GetVolume(), &hmat1, &hmat2, kTRUE, ovlp);
1825
1826 next2.Skip();
1827 }
1828 }
1829
1830 } else {
1831 //=== ASSEMBLY/REAL <-> REAL
1832 // Selected node skip logic
1833 if (fSelectedNode) {
1834 if ((fSelectedNode != node) && (fSelectedNode != node02) &&
1836 next1.Skip();
1837 continue;
1838 }
1839 if ((fSelectedNode != node) && (fSelectedNode != node02)) {
1840 level = next1.GetLevel();
1841 while ((nodecheck = next1.GetNode(level--))) {
1842 if (nodecheck == fSelectedNode)
1843 break;
1844 }
1845 if (node01 == fSelectedNode)
1846 nodecheck = node01;
1847 if (!nodecheck) {
1848 next1.Skip();
1849 continue;
1850 }
1851 }
1852 }
1853
1854 ncand++;
1855 PushCandidate(out,
1856 TString::Format("%s/%s overlapping %s/%s", vol->GetName(), path.Data(),
1857 vol->GetName(), path1.Data()),
1858 node->GetVolume(), node02->GetVolume(), &hmat1, node02->GetMatrix(), kTRUE, ovlp);
1859 }
1860
1861 next1.Skip();
1862 }
1863 }
1864
1865 } else {
1866 if (node02->GetVolume()->IsAssembly()) {
1867 next2.Reset(node02->GetVolume());
1868 //=== REAL <-> ASSEMBLY/ANY ===//
1869 while ((node1 = next2())) {
1870 hmat2 = node02->GetMatrix();
1871 hmat2 *= *next2.GetCurrentMatrix();
1872 box2 = (const TGeoBBox *)node1->GetVolume()->GetShape();
1873 // OBB check
1874 if (!TGeoBBox::MayIntersect(box1, node01->GetMatrix(), box2, &hmat2)) {
1875 // No intersection, skip the entire right branch
1876 next2.Skip();
1877 continue;
1878 }
1879
1880 if (!node1->GetVolume()->IsAssembly()) {
1881 // Selected node skip logic
1882 if (fSelectedNode) {
1883 if ((fSelectedNode != node1) && (fSelectedNode != node01) &&
1885 next2.Skip();
1886 continue;
1887 }
1888 if ((fSelectedNode != node1) && (fSelectedNode != node01)) {
1889 level = next2.GetLevel();
1890 while ((nodecheck = next2.GetNode(level--))) {
1891 if (nodecheck == fSelectedNode)
1892 break;
1893 }
1894 if (node02 == fSelectedNode)
1895 nodecheck = node02;
1896 if (!nodecheck) {
1897 next2.Skip();
1898 continue;
1899 }
1900 }
1901 }
1902
1903 ncand++;
1904 next2.GetPath(path1);
1905 PushCandidate(out,
1906 TString::Format("%s/%s overlapping %s/%s", vol->GetName(), path.Data(),
1907 vol->GetName(), path1.Data()),
1908 node01->GetVolume(), node1->GetVolume(), node01->GetMatrix(), &hmat2, kTRUE, ovlp);
1909
1910 next2.Skip();
1911 }
1912 }
1913
1914 } else {
1915 // both not assembly
1917 continue;
1918
1919 ncand++;
1921 out,
1922 TString::Format("%s/%s overlapping %s/%s", vol->GetName(), path.Data(), vol->GetName(), path1.Data()),
1923 node01->GetVolume(), node02->GetVolume(), node01->GetMatrix(), node02->GetMatrix(), kTRUE, ovlp);
1924 }
1925 }
1926 }
1927 }
1928 // if ((ncand - nextr) > 0)
1929 // printf("overlap cand for %s : %d\n", vol->GetName(), ncand - nextr);
1930 return ncand;
1931}
1932
1933////////////////////////////////////////////////////////////////////////////////
1934/// Helper to fill candidates list
1935
1936void TGeoChecker::PushCandidate(std::vector<TGeoOverlapCandidate> &out, const TString &name, TGeoVolume *vol1,
1938 Double_t ovlp) const
1939{
1941 c.fName = name;
1942 c.fVol1 = vol1;
1943 c.fVol2 = vol2;
1944 c.fMat1 = TGeoHMatrix(*mat1);
1945 c.fMat2 = TGeoHMatrix(*mat2);
1946 c.fIsOverlap = isovlp;
1947 c.fOvlp = ovlp;
1948 out.emplace_back(std::move(c));
1949}
1950
1951////////////////////////////////////////////////////////////////////////////////
1952/// Print the current list of overlaps held by the manager class.
1953
1955{
1957 TGeoOverlap *ov;
1958 printf("=== Overlaps for %s ===\n", fGeoManager->GetName());
1959 while ((ov = (TGeoOverlap *)next()))
1960 ov->PrintInfo();
1961}
1962
1963////////////////////////////////////////////////////////////////////////////////
1964/// Draw point (x,y,z) over the picture of the daughters of the volume containing this point.
1965/// Generates a report regarding the path to the node containing this point and the distance to
1966/// the closest boundary.
1967
1969{
1970 Double_t point[3];
1971 Double_t local[3];
1972 point[0] = x;
1973 point[1] = y;
1974 point[2] = z;
1976 if (fVsafe) {
1977 TGeoNode *old = fVsafe->GetNode("SAFETY_1");
1978 if (old)
1979 fVsafe->GetNodes()->RemoveAt(vol->GetNdaughters() - 1);
1980 }
1981 // if (vol != fGeoManager->GetMasterVolume()) fGeoManager->RestoreMasterVolume();
1982 TGeoNode *node = fGeoManager->FindNode(point[0], point[1], point[2]);
1984 // get current node
1985 printf("=== Check current point : (%g, %g, %g) ===\n", point[0], point[1], point[2]);
1986 printf(" - path : %s\n", fGeoManager->GetPath());
1987 // get corresponding volume
1988 if (node)
1989 vol = node->GetVolume();
1990 // compute safety distance (distance to boundary ignored)
1991 Double_t close = (safety > 0.) ? safety : fGeoManager->Safety();
1992 printf("Safety radius : %f\n", close);
1993 if (close > 1E-4) {
1994 TGeoVolume *sph = fGeoManager->MakeSphere("SAFETY", vol->GetMedium(), 0, close, 0, 180, 0, 360);
1995 sph->SetLineColor(2);
1996 sph->SetLineStyle(3);
1997 vol->AddNode(sph, 1, new TGeoTranslation(local[0], local[1], local[2]));
1998 fVsafe = vol;
1999 }
2001 pm->SetMarkerColor(2);
2002 pm->SetMarkerStyle(8);
2003 pm->SetMarkerSize(0.5);
2004 pm->SetNextPoint(local[0], local[1], local[2]);
2005 if (vol->GetNdaughters() < 2)
2007 else
2010 if (!vol->IsVisible())
2011 vol->SetVisibility(kTRUE);
2012 vol->Draw();
2013 pm->Draw("SAME");
2014 gPad->Modified();
2015 gPad->Update();
2016}
2017
2018////////////////////////////////////////////////////////////////////////////////
2019/// Test for shape navigation methods. Summary for test numbers:
2020/// - 1: DistFromInside/Outside. Sample points inside the shape. Generate
2021/// directions randomly in cos(theta). Compute DistFromInside and move the
2022/// point with bigger distance. Compute DistFromOutside back from new point.
2023/// Plot d-(d1+d2)
2024/// - 2: Safety test. Sample points inside the bounding and compute safety. Generate
2025/// directions randomly in cos(theta) and compute distance to boundary. Check if
2026/// distance to boundary is bigger than safety
2027
2029{
2030 switch (testNo) {
2031 case 1: ShapeDistances(shape, nsamples, option); break;
2032 case 2: ShapeSafety(shape, nsamples, option); break;
2033 case 3: ShapeNormal(shape, nsamples, option); break;
2034 default: Error("CheckShape", "Test number %d not existent", testNo);
2035 }
2036}
2037
2038////////////////////////////////////////////////////////////////////////////////
2039/// Test TGeoShape::DistFromInside/Outside. Sample points inside the shape. Generate
2040/// directions randomly in cos(theta). Compute d1 = DistFromInside and move the
2041/// point on the boundary. Compute DistFromOutside and propagate with d2 making sure that
2042/// the shape is not re-entered. Swap direction and call DistFromOutside that
2043/// should fall back on the same point on the boundary (at d2). Propagate back on boundary
2044/// then compute DistFromInside that should be bigger than d1.
2045/// Plot d-(d1+d2)
2046
2048{
2049 Double_t dx = ((TGeoBBox *)shape)->GetDX();
2050 Double_t dy = ((TGeoBBox *)shape)->GetDY();
2051 Double_t dz = ((TGeoBBox *)shape)->GetDZ();
2052 Double_t dmax = 2. * TMath::Sqrt(dx * dx + dy * dy + dz * dz);
2054 Int_t itot = 0;
2055 // Number of tracks shot for every point inside the shape
2056 const Int_t kNtracks = 1000;
2057 Int_t n10 = nsamples / 10;
2058 Int_t i, j;
2059 Double_t point[3], pnew[3];
2060 Double_t dir[3], dnew[3];
2061 Double_t theta, phi, delta;
2062 TPolyMarker3D *pmfrominside = nullptr;
2063 TPolyMarker3D *pmfromoutside = nullptr;
2064 TH1D *hist = new TH1D("hTest1", "Residual distance from inside/outside", 200, -20, 0);
2065 hist->GetXaxis()->SetTitle("delta[cm] - first bin=overflow");
2066 hist->GetYaxis()->SetTitle("count");
2068
2069 if (!fTimer)
2070 fTimer = new TStopwatch();
2071 fTimer->Reset();
2072 fTimer->Start();
2073 while (itot < nsamples) {
2074 Bool_t inside = kFALSE;
2075 while (!inside) {
2076 point[0] = gRandom->Uniform(-dx, dx);
2077 point[1] = gRandom->Uniform(-dy, dy);
2078 point[2] = gRandom->Uniform(-dz, dz);
2079 inside = shape->Contains(point);
2080 }
2081 itot++;
2082 if (n10) {
2083 if ((itot % n10) == 0)
2084 printf("%i percent\n", Int_t(100 * itot / nsamples));
2085 }
2086 for (i = 0; i < kNtracks; i++) {
2087 phi = 2 * TMath::Pi() * gRandom->Rndm();
2088 theta = TMath::ACos(1. - 2. * gRandom->Rndm());
2089 dir[0] = TMath::Sin(theta) * TMath::Cos(phi);
2090 dir[1] = TMath::Sin(theta) * TMath::Sin(phi);
2091 dir[2] = TMath::Cos(theta);
2092 dmove = dmax;
2093 // We have track direction, compute distance from inside
2094 d1 = shape->DistFromInside(point, dir, 3);
2095 if (d1 > dmove || d1 < TGeoShape::Tolerance()) {
2096 // Bad distance or bbox size, to debug
2097 new TCanvas("shape01", Form("Shape %s (%s)", shape->GetName(), shape->ClassName()), 1000, 800);
2098 shape->Draw();
2099 printf("DistFromInside: (%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) %f/%f(max)\n", point[0],
2100 point[1], point[2], dir[0], dir[1], dir[2], d1, dmove);
2101 pmfrominside = new TPolyMarker3D(2);
2102 pmfrominside->SetMarkerColor(kRed);
2103 pmfrominside->SetMarkerStyle(24);
2104 pmfrominside->SetMarkerSize(0.4);
2105 pmfrominside->SetNextPoint(point[0], point[1], point[2]);
2106 for (j = 0; j < 3; j++)
2107 pnew[j] = point[j] + d1 * dir[j];
2108 pmfrominside->SetNextPoint(pnew[0], pnew[1], pnew[2]);
2109 pmfrominside->Draw();
2110 return;
2111 }
2112 // Propagate BEFORE the boundary and make sure that DistFromOutside
2113 // does not return 0 (!!!)
2114 // Check if there is a second crossing
2115 for (j = 0; j < 3; j++)
2116 pnew[j] = point[j] + (d1 - TGeoShape::Tolerance()) * dir[j];
2117 dnext = shape->DistFromOutside(pnew, dir, 3);
2118 if (d1 + dnext < dmax)
2119 dmove = d1 + 0.5 * dnext;
2120 // Move point and swap direction
2121 for (j = 0; j < 3; j++) {
2122 pnew[j] = point[j] + dmove * dir[j];
2123 dnew[j] = -dir[j];
2124 }
2125 // Compute now distance from outside
2126 d2 = shape->DistFromOutside(pnew, dnew, 3);
2127 delta = dmove - d1 - d2;
2128 if (TMath::Abs(delta) > 1E-6 || dnext < 2. * TGeoShape::Tolerance()) {
2129 // Error->debug this
2130 new TCanvas("shape01", Form("Shape %s (%s)", shape->GetName(), shape->ClassName()), 1000, 800);
2131 shape->Draw();
2132 printf("Error: (%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) d1=%f d2=%f dmove=%f\n", point[0],
2133 point[1], point[2], dir[0], dir[1], dir[2], d1, d2, dmove);
2134 if (dnext < 2. * TGeoShape::Tolerance()) {
2135 printf(" (*)DistFromOutside(%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) dnext = %f\n",
2136 point[0] + (d1 - TGeoShape::Tolerance()) * dir[0],
2137 point[1] + (d1 - TGeoShape::Tolerance()) * dir[1],
2138 point[2] + (d1 - TGeoShape::Tolerance()) * dir[2], dir[0], dir[1], dir[2], dnext);
2139 } else {
2140 printf(" DistFromOutside(%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) dnext = %f\n",
2141 point[0] + d1 * dir[0], point[1] + d1 * dir[1], point[2] + d1 * dir[2], dir[0], dir[1], dir[2],
2142 dnext);
2143 }
2144 printf(" DistFromOutside(%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) = %f\n", pnew[0], pnew[1],
2145 pnew[2], dnew[0], dnew[1], dnew[2], d2);
2146 pmfrominside = new TPolyMarker3D(2);
2147 pmfrominside->SetMarkerStyle(24);
2148 pmfrominside->SetMarkerSize(0.4);
2149 pmfrominside->SetMarkerColor(kRed);
2150 pmfrominside->SetNextPoint(point[0], point[1], point[2]);
2151 for (j = 0; j < 3; j++)
2152 point[j] += d1 * dir[j];
2153 pmfrominside->SetNextPoint(point[0], point[1], point[2]);
2154 pmfrominside->Draw();
2156 pmfromoutside->SetMarkerStyle(20);
2157 pmfromoutside->SetMarkerStyle(7);
2158 pmfromoutside->SetMarkerSize(0.3);
2159 pmfromoutside->SetMarkerColor(kBlue);
2160 pmfromoutside->SetNextPoint(pnew[0], pnew[1], pnew[2]);
2161 for (j = 0; j < 3; j++)
2162 pnew[j] += d2 * dnew[j];
2163 if (d2 < 1E10)
2164 pmfromoutside->SetNextPoint(pnew[0], pnew[1], pnew[2]);
2165 pmfromoutside->Draw();
2166 return;
2167 }
2168 // Compute distance from inside which should be bigger than d1
2169 for (j = 0; j < 3; j++)
2170 pnew[j] += (d2 - TGeoShape::Tolerance()) * dnew[j];
2171 dnext = shape->DistFromInside(pnew, dnew, 3);
2172 if (dnext < d1 - TGeoShape::Tolerance() || dnext > dmax) {
2173 new TCanvas("shape01", Form("Shape %s (%s)", shape->GetName(), shape->ClassName()), 1000, 800);
2174 shape->Draw();
2175 printf("Error DistFromInside(%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) d1=%f d1p=%f\n", pnew[0],
2176 pnew[1], pnew[2], dnew[0], dnew[1], dnew[2], d1, dnext);
2177 pmfrominside = new TPolyMarker3D(2);
2178 pmfrominside->SetMarkerStyle(24);
2179 pmfrominside->SetMarkerSize(0.4);
2180 pmfrominside->SetMarkerColor(kRed);
2181 pmfrominside->SetNextPoint(point[0], point[1], point[2]);
2182 for (j = 0; j < 3; j++)
2183 point[j] += d1 * dir[j];
2184 pmfrominside->SetNextPoint(point[0], point[1], point[2]);
2185 pmfrominside->Draw();
2187 pmfromoutside->SetMarkerStyle(20);
2188 pmfromoutside->SetMarkerStyle(7);
2189 pmfromoutside->SetMarkerSize(0.3);
2190 pmfromoutside->SetMarkerColor(kBlue);
2191 pmfromoutside->SetNextPoint(pnew[0], pnew[1], pnew[2]);
2192 for (j = 0; j < 3; j++)
2193 pnew[j] += dnext * dnew[j];
2194 if (d2 < 1E10)
2195 pmfromoutside->SetNextPoint(pnew[0], pnew[1], pnew[2]);
2196 pmfromoutside->Draw();
2197 return;
2198 }
2199 if (TMath::Abs(delta) < 1E-20)
2200 delta = 1E-30;
2201 hist->Fill(TMath::Max(TMath::Log(TMath::Abs(delta)), -20.));
2202 }
2203 }
2204 fTimer->Stop();
2205 fTimer->Print();
2206 new TCanvas("Test01", "Residuals DistFromInside/Outside", 800, 600);
2207 hist->Draw();
2208}
2209
2210////////////////////////////////////////////////////////////////////////////////
2211/// Check of validity of safe distance for a given shape.
2212/// Sample points inside the 2x bounding box and compute safety. Generate
2213/// directions randomly in cos(theta) and compute distance to boundary. Check if
2214/// distance to boundary is bigger than safety.
2215
2217{
2218 Double_t dx = ((TGeoBBox *)shape)->GetDX();
2219 Double_t dy = ((TGeoBBox *)shape)->GetDY();
2220 Double_t dz = ((TGeoBBox *)shape)->GetDZ();
2221 // Number of tracks shot for every point inside the shape
2222 const Int_t kNtracks = 1000;
2223 Int_t n10 = nsamples / 10;
2224 Int_t i;
2225 Double_t dist;
2226 Double_t point[3];
2227 Double_t dir[3];
2228 Double_t theta, phi;
2229 TPolyMarker3D *pm1 = nullptr;
2230 TPolyMarker3D *pm2 = nullptr;
2231 if (!fTimer)
2232 fTimer = new TStopwatch();
2233 fTimer->Reset();
2234 fTimer->Start();
2235 Int_t itot = 0;
2236 while (itot < nsamples) {
2237 Bool_t inside = kFALSE;
2238 point[0] = gRandom->Uniform(-2 * dx, 2 * dx);
2239 point[1] = gRandom->Uniform(-2 * dy, 2 * dy);
2240 point[2] = gRandom->Uniform(-2 * dz, 2 * dz);
2241 inside = shape->Contains(point);
2242 Double_t safe = shape->Safety(point, inside);
2243 itot++;
2244 if (n10) {
2245 if ((itot % n10) == 0)
2246 printf("%i percent\n", Int_t(100 * itot / nsamples));
2247 }
2248 for (i = 0; i < kNtracks; i++) {
2249 phi = 2 * TMath::Pi() * gRandom->Rndm();
2250 theta = TMath::ACos(1. - 2. * gRandom->Rndm());
2251 dir[0] = TMath::Sin(theta) * TMath::Cos(phi);
2252 dir[1] = TMath::Sin(theta) * TMath::Sin(phi);
2253 dir[2] = TMath::Cos(theta);
2254 if (inside)
2255 dist = shape->DistFromInside(point, dir, 3);
2256 else
2257 dist = shape->DistFromOutside(point, dir, 3);
2258 if (dist < safe) {
2259 printf("Error safety (%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) safe=%f dist=%f\n", point[0],
2260 point[1], point[2], dir[0], dir[1], dir[2], safe, dist);
2261 new TCanvas("shape02", Form("Shape %s (%s)", shape->GetName(), shape->ClassName()), 1000, 800);
2262 shape->Draw();
2263 pm1 = new TPolyMarker3D(2);
2264 pm1->SetMarkerStyle(24);
2265 pm1->SetMarkerSize(0.4);
2266 pm1->SetMarkerColor(kRed);
2267 pm1->SetNextPoint(point[0], point[1], point[2]);
2268 pm1->SetNextPoint(point[0] + safe * dir[0], point[1] + safe * dir[1], point[2] + safe * dir[2]);
2269 pm1->Draw();
2270 pm2 = new TPolyMarker3D(1);
2271 pm2->SetMarkerStyle(7);
2272 pm2->SetMarkerSize(0.3);
2273 pm2->SetMarkerColor(kBlue);
2274 pm2->SetNextPoint(point[0] + dist * dir[0], point[1] + dist * dir[1], point[2] + dist * dir[2]);
2275 pm2->Draw();
2276 return;
2277 }
2278 }
2279 }
2280 fTimer->Stop();
2281 fTimer->Print();
2282}
2283
2284////////////////////////////////////////////////////////////////////////////////
2285/// Check of validity of the normal for a given shape.
2286/// Sample points inside the shape. Generate directions randomly in cos(theta)
2287/// and propagate to boundary. Compute normal and safety at crossing point, plot
2288/// the point and generate a random direction so that (dir) dot (norm) <0.
2289
2291{
2292 Double_t dx = ((TGeoBBox *)shape)->GetDX();
2293 Double_t dy = ((TGeoBBox *)shape)->GetDY();
2294 Double_t dz = ((TGeoBBox *)shape)->GetDZ();
2295 Double_t dmax = 2. * TMath::Sqrt(dx * dx + dy * dy + dz * dz);
2296 // Number of tracks shot for every point inside the shape
2297 const Int_t kNtracks = 1000;
2298 Int_t n10 = nsamples / 10;
2299 Int_t itot = 0, errcnt = 0, errsame = 0;
2300 Int_t i;
2301 Double_t dist, olddist, safe, dot;
2302 Double_t theta, phi, ndotd;
2303 Double_t *spoint = new Double_t[3 * nsamples];
2304 Double_t *sdir = new Double_t[3 * nsamples];
2305 while (itot < nsamples) {
2306 Bool_t inside = kFALSE;
2307 while (!inside) {
2308 spoint[3 * itot] = gRandom->Uniform(-dx, dx);
2309 spoint[3 * itot + 1] = gRandom->Uniform(-dy, dy);
2310 spoint[3 * itot + 2] = gRandom->Uniform(-dz, dz);
2311 inside = shape->Contains(&spoint[3 * itot]);
2312 }
2313 phi = 2 * TMath::Pi() * gRandom->Rndm();
2314 theta = TMath::ACos(1. - 2. * gRandom->Rndm());
2315 sdir[3 * itot] = TMath::Sin(theta) * TMath::Cos(phi);
2316 sdir[3 * itot + 1] = TMath::Sin(theta) * TMath::Sin(phi);
2317 sdir[3 * itot + 2] = TMath::Cos(theta);
2318 itot++;
2319 }
2320 Double_t point[3], newpoint[3], oldpoint[3];
2321 Double_t dir[3], olddir[3];
2322 Double_t norm[3], newnorm[3], oldnorm[3];
2323 TCanvas *errcanvas = nullptr;
2324 TPolyMarker3D *pm1 = nullptr;
2325 TPolyMarker3D *pm2 = nullptr;
2326 pm2 = new TPolyMarker3D();
2327 // pm2->SetMarkerStyle(24);
2328 pm2->SetMarkerSize(0.2);
2329 pm2->SetMarkerColor(kBlue);
2330 if (!fTimer)
2331 fTimer = new TStopwatch();
2332 fTimer->Reset();
2333 fTimer->Start();
2334 for (itot = 0; itot < nsamples; itot++) {
2335 if (n10) {
2336 if ((itot % n10) == 0)
2337 printf("%i percent\n", Int_t(100 * itot / nsamples));
2338 }
2339 oldnorm[0] = oldnorm[1] = oldnorm[2] = 0.;
2340 olddist = 0.;
2341 for (Int_t j = 0; j < 3; j++) {
2342 oldpoint[j] = point[j] = spoint[3 * itot + j];
2343 olddir[j] = dir[j] = sdir[3 * itot + j];
2344 }
2345 for (i = 0; i < kNtracks; i++) {
2346 if (errcnt > 0)
2347 break;
2348 dist = shape->DistFromInside(point, dir, 3);
2349 for (Int_t j = 0; j < 3; j++) {
2350 newpoint[j] = point[j] + dist * dir[j];
2351 }
2352 shape->ComputeNormal(newpoint, dir, newnorm);
2353
2354 dot = olddir[0] * oldnorm[0] + olddir[1] * oldnorm[1] + olddir[2] * oldnorm[2];
2355 if (!shape->Contains(point) && shape->Safety(point, kFALSE) > 1.E-3) {
2356 errcnt++;
2357 printf("Error point outside (%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) =%g olddist=%g\n",
2358 point[0], point[1], point[2], dir[0], dir[1], dir[2], dist, olddist);
2359 printf(" old point: (%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f)\n", oldpoint[0],
2360 oldpoint[1], oldpoint[2], olddir[0], olddir[1], olddir[2]);
2361 if (!errcanvas)
2362 errcanvas =
2363 new TCanvas("shape_err03", Form("Shape %s (%s)", shape->GetName(), shape->ClassName()), 1000, 800);
2364 if (!pm1) {
2365 pm1 = new TPolyMarker3D();
2366 pm1->SetMarkerStyle(24);
2367 pm1->SetMarkerSize(0.4);
2368 pm1->SetMarkerColor(kRed);
2369 }
2370 pm1->SetNextPoint(point[0], point[1], point[2]);
2371 pm1->SetNextPoint(oldpoint[0], oldpoint[1], oldpoint[2]);
2372 break;
2373 }
2374 if ((dist < TGeoShape::Tolerance() && olddist * dot > 1.E-3) || dist > dmax) {
2375 errsame++;
2376 if (errsame > 1) {
2377 errcnt++;
2378 printf("Error DistFromInside(%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f) =%g olddist=%g\n",
2379 point[0], point[1], point[2], dir[0], dir[1], dir[2], dist, olddist);
2380 printf(" new norm: (%g, %g, %g)\n", newnorm[0], newnorm[1], newnorm[2]);
2381 printf(" old point: (%19.15f, %19.15f, %19.15f, %19.15f, %19.15f, %19.15f)\n", oldpoint[0],
2382 oldpoint[1], oldpoint[2], olddir[0], olddir[1], olddir[2]);
2383 printf(" old norm: (%g, %g, %g)\n", oldnorm[0], oldnorm[1], oldnorm[2]);
2384 if (!errcanvas)
2385 errcanvas =
2386 new TCanvas("shape_err03", Form("Shape %s (%s)", shape->GetName(), shape->ClassName()), 1000, 800);
2387 if (!pm1) {
2388 pm1 = new TPolyMarker3D();
2389 pm1->SetMarkerStyle(24);
2390 pm1->SetMarkerSize(0.4);
2391 pm1->SetMarkerColor(kRed);
2392 }
2393 pm1->SetNextPoint(point[0], point[1], point[2]);
2394 pm1->SetNextPoint(oldpoint[0], oldpoint[1], oldpoint[2]);
2395 break;
2396 }
2397 } else
2398 errsame = 0;
2399
2400 olddist = dist;
2401 for (Int_t j = 0; j < 3; j++) {
2402 oldpoint[j] = point[j];
2403 point[j] += dist * dir[j];
2404 }
2405 safe = shape->Safety(point, kTRUE);
2406 if (safe > 1.E-3) {
2407 errcnt++;
2408 printf("Error safety (%19.15f, %19.15f, %19.15f) safe=%g\n", point[0], point[1], point[2], safe);
2409 if (!errcanvas)
2410 errcanvas =
2411 new TCanvas("shape_err03", Form("Shape %s (%s)", shape->GetName(), shape->ClassName()), 1000, 800);
2412 if (!pm1) {
2413 pm1 = new TPolyMarker3D();
2414 pm1->SetMarkerStyle(24);
2415 pm1->SetMarkerSize(0.4);
2416 pm1->SetMarkerColor(kRed);
2417 }
2418 pm1->SetNextPoint(point[0], point[1], point[2]);
2419 break;
2420 }
2421 // Compute normal
2422 shape->ComputeNormal(point, dir, norm);
2423 memcpy(oldnorm, norm, 3 * sizeof(Double_t));
2424 memcpy(olddir, dir, 3 * sizeof(Double_t));
2425 while (true) {
2426 phi = 2 * TMath::Pi() * gRandom->Rndm();
2427 theta = TMath::ACos(1. - 2. * gRandom->Rndm());
2428 dir[0] = TMath::Sin(theta) * TMath::Cos(phi);
2429 dir[1] = TMath::Sin(theta) * TMath::Sin(phi);
2430 dir[2] = TMath::Cos(theta);
2431 ndotd = dir[0] * norm[0] + dir[1] * norm[1] + dir[2] * norm[2];
2432 if (ndotd < 0)
2433 break; // backwards, still inside shape
2434 }
2435 if ((itot % 10) == 0)
2436 pm2->SetNextPoint(point[0], point[1], point[2]);
2437 }
2438 }
2439 fTimer->Stop();
2440 fTimer->Print();
2441 if (errcanvas) {
2442 shape->Draw();
2443 pm1->Draw();
2444 }
2445
2446 new TCanvas("shape03", Form("Shape %s (%s)", shape->GetName(), shape->ClassName()), 1000, 800);
2447 pm2->Draw();
2448 delete[] spoint;
2449 delete[] sdir;
2450}
2451
2452////////////////////////////////////////////////////////////////////////////////
2453/// Generate a lego plot fot the top volume, according to option.
2454
2456 Double_t phimax, Double_t /*rmin*/, Double_t /*rmax*/, Option_t *option)
2457{
2458 TH2F *hist = new TH2F("lego", option, nphi, phimin, phimax, ntheta, themin, themax);
2459
2460 Double_t degrad = TMath::Pi() / 180.;
2461 Double_t theta, phi, step, matprop, x;
2462 Double_t start[3];
2463 Double_t dir[3];
2465 Int_t i; // loop index for phi
2466 Int_t j; // loop index for theta
2467 Int_t ntot = ntheta * nphi;
2468 Int_t n10 = ntot / 10;
2469 Int_t igen = 0, iloop = 0;
2470 printf("=== Lego plot sph. => nrays=%i\n", ntot);
2471 for (i = 1; i <= nphi; i++) {
2472 for (j = 1; j <= ntheta; j++) {
2473 igen++;
2474 if (n10) {
2475 if ((igen % n10) == 0)
2476 printf("%i percent\n", Int_t(100 * igen / ntot));
2477 }
2478 x = 0;
2479 theta = hist->GetYaxis()->GetBinCenter(j);
2480 phi = hist->GetXaxis()->GetBinCenter(i) + 1E-3;
2481 start[0] = start[1] = start[2] = 1E-3;
2482 dir[0] = TMath::Sin(theta * degrad) * TMath::Cos(phi * degrad);
2483 dir[1] = TMath::Sin(theta * degrad) * TMath::Sin(phi * degrad);
2484 dir[2] = TMath::Cos(theta * degrad);
2485 fGeoManager->InitTrack(&start[0], &dir[0]);
2487 if (fGeoManager->IsOutside())
2488 startnode = nullptr;
2489 if (startnode) {
2490 matprop = startnode->GetVolume()->GetMaterial()->GetRadLen();
2491 } else {
2492 matprop = 0.;
2493 }
2495 // fGeoManager->IsStepEntering();
2496 // find where we end-up
2498 step = fGeoManager->GetStep();
2499 while (step < 1E10) {
2500 // now see if we can make an other step
2501 iloop = 0;
2502 while (!fGeoManager->IsEntering()) {
2503 iloop++;
2504 fGeoManager->SetStep(1E-3);
2505 step += 1E-3;
2507 }
2508 if (iloop > 1000)
2509 printf("%i steps\n", iloop);
2510 if (matprop > 0) {
2511 x += step / matprop;
2512 }
2513 if (endnode == nullptr && step > 1E10)
2514 break;
2515 // generate an extra step to cross boundary
2517 if (startnode) {
2518 matprop = startnode->GetVolume()->GetMaterial()->GetRadLen();
2519 } else {
2520 matprop = 0.;
2521 }
2522
2525 step = fGeoManager->GetStep();
2526 }
2527 hist->Fill(phi, theta, x);
2528 }
2529 }
2530 return hist;
2531}
2532
2533////////////////////////////////////////////////////////////////////////////////
2534/// Draw random points in the bounding box of a volume.
2535
2537{
2538 if (!vol)
2539 return;
2540 vol->VisibleDaughters(kTRUE);
2541 vol->Draw();
2542 TString opt = option;
2543 opt.ToLower();
2544 TObjArray *pm = new TObjArray(128);
2545 TPolyMarker3D *marker = nullptr;
2546 const TGeoShape *shape = vol->GetShape();
2547 TGeoBBox *box = (TGeoBBox *)shape;
2548 Double_t dx = box->GetDX();
2549 Double_t dy = box->GetDY();
2550 Double_t dz = box->GetDZ();
2551 Double_t ox = (box->GetOrigin())[0];
2552 Double_t oy = (box->GetOrigin())[1];
2553 Double_t oz = (box->GetOrigin())[2];
2554 Double_t *xyz = new Double_t[3];
2555 printf("Random box : %f, %f, %f\n", dx, dy, dz);
2556 TGeoNode *node = nullptr;
2557 printf("Start... %i points\n", npoints);
2558 Int_t i = 0;
2559 Int_t igen = 0;
2560 Int_t ic = 0;
2561 Int_t n10 = npoints / 10;
2562 Double_t ratio = 0;
2563 while (igen < npoints) {
2564 xyz[0] = ox - dx + 2 * dx * gRandom->Rndm();
2565 xyz[1] = oy - dy + 2 * dy * gRandom->Rndm();
2566 xyz[2] = oz - dz + 2 * dz * gRandom->Rndm();
2568 igen++;
2569 if (n10) {
2570 if ((igen % n10) == 0)
2571 printf("%i percent\n", Int_t(100 * igen / npoints));
2572 }
2573 node = fGeoManager->FindNode();
2574 if (!node)
2575 continue;
2576 if (!node->IsOnScreen())
2577 continue;
2578 // draw only points in overlapping/non-overlapping volumes
2579 if (opt.Contains("many") && !node->IsOverlapping())
2580 continue;
2581 if (opt.Contains("only") && node->IsOverlapping())
2582 continue;
2583 ic = node->GetColour();
2584 if ((ic < 0) || (ic >= 128))
2585 ic = 1;
2586 marker = (TPolyMarker3D *)pm->At(ic);
2587 if (!marker) {
2588 marker = new TPolyMarker3D();
2589 marker->SetMarkerColor(ic);
2590 // marker->SetMarkerStyle(8);
2591 // marker->SetMarkerSize(0.4);
2592 pm->AddAt(marker, ic);
2593 }
2594 marker->SetNextPoint(xyz[0], xyz[1], xyz[2]);
2595 i++;
2596 }
2597 printf("Number of visible points : %i\n", i);
2598 if (igen)
2599 ratio = (Double_t)i / (Double_t)igen;
2600 printf("efficiency : %g\n", ratio);
2601 for (Int_t m = 0; m < 128; m++) {
2602 marker = (TPolyMarker3D *)pm->At(m);
2603 if (marker)
2604 marker->Draw("SAME");
2605 }
2607 printf("---Daughters of %s made invisible.\n", fGeoManager->GetTopVolume()->GetName());
2608 printf("---Make them visible with : gGeoManager->GetTopVolume()->VisibleDaughters();\n");
2609 delete pm;
2610 delete[] xyz;
2611}
2612
2613////////////////////////////////////////////////////////////////////////////////
2614/// Randomly shoot nrays from point (startx,starty,startz) and plot intersections
2615/// with surfaces for current top node.
2616
2619{
2620 TObjArray *pm = new TObjArray(128);
2622 TPolyLine3D *line = nullptr;
2623 TPolyLine3D *normline = nullptr;
2625 // vol->VisibleDaughters(kTRUE);
2626
2628 if (nrays <= 0) {
2629 nrays = 100000;
2630 random = kTRUE;
2631 }
2632 Double_t start[3];
2633 Double_t dir[3];
2634 const Double_t *point = fGeoManager->GetCurrentPoint();
2635 vol->Draw();
2636 printf("Start... %i rays\n", nrays);
2638 Bool_t vis1, vis2;
2639 Int_t i = 0;
2641 Int_t itot = 0;
2642 Int_t n10 = nrays / 10;
2643 Double_t theta, phi, step, normlen;
2644 Double_t ox = ((TGeoBBox *)vol->GetShape())->GetOrigin()[0];
2645 Double_t oy = ((TGeoBBox *)vol->GetShape())->GetOrigin()[1];
2646 Double_t oz = ((TGeoBBox *)vol->GetShape())->GetOrigin()[2];
2647 Double_t dx = ((TGeoBBox *)vol->GetShape())->GetDX();
2648 Double_t dy = ((TGeoBBox *)vol->GetShape())->GetDY();
2649 Double_t dz = ((TGeoBBox *)vol->GetShape())->GetDZ();
2650 normlen = TMath::Max(dx, dy);
2652 normlen *= 0.05;
2653 while (itot < nrays) {
2654 itot++;
2655 inull = 0;
2656 ipoint = 0;
2657 if (n10) {
2658 if ((itot % n10) == 0)
2659 printf("%i percent\n", Int_t(100 * itot / nrays));
2660 }
2661 if (random) {
2662 start[0] = ox - dx + 2 * dx * gRandom->Rndm();
2663 start[1] = oy - dy + 2 * dy * gRandom->Rndm();
2664 start[2] = oz - dz + 2 * dz * gRandom->Rndm();
2665 } else {
2666 start[0] = startx;
2667 start[1] = starty;
2668 start[2] = startz;
2669 }
2670 phi = 2 * TMath::Pi() * gRandom->Rndm();
2671 theta = TMath::ACos(1. - 2. * gRandom->Rndm());
2672 dir[0] = TMath::Sin(theta) * TMath::Cos(phi);
2673 dir[1] = TMath::Sin(theta) * TMath::Sin(phi);
2674 dir[2] = TMath::Cos(theta);
2675 startnode = fGeoManager->InitTrack(start[0], start[1], start[2], dir[0], dir[1], dir[2]);
2676 line = nullptr;
2677 if (fGeoManager->IsOutside())
2678 startnode = nullptr;
2679 vis1 = kFALSE;
2680 if (target_vol) {
2681 if (startnode && starget == startnode->GetVolume()->GetName())
2682 vis1 = kTRUE;
2683 } else {
2684 if (startnode && startnode->IsOnScreen())
2685 vis1 = kTRUE;
2686 }
2687 if (vis1) {
2688 line = new TPolyLine3D(2);
2689 line->SetLineColor(startnode->GetVolume()->GetLineColor());
2690 line->SetPoint(ipoint++, start[0], start[1], start[2]);
2691 i++;
2692 pm->Add(line);
2693 }
2695 step = fGeoManager->GetStep();
2696 if (step < TGeoShape::Tolerance())
2697 inull++;
2698 else
2699 inull = 0;
2700 if (inull > 5)
2701 break;
2702 const Double_t *normal = nullptr;
2703 if (check_norm) {
2705 if (!normal)
2706 break;
2707 }
2708 vis2 = kFALSE;
2709 if (target_vol) {
2710 if (starget == endnode->GetVolume()->GetName())
2711 vis2 = kTRUE;
2712 } else if (endnode->IsOnScreen())
2713 vis2 = kTRUE;
2714 if (ipoint > 0) {
2715 // old visible node had an entry point -> finish segment
2716 line->SetPoint(ipoint, point[0], point[1], point[2]);
2717 if (!vis2 && check_norm) {
2718 normline = new TPolyLine3D(2);
2719 normline->SetLineColor(kBlue);
2720 normline->SetLineWidth(1);
2721 normline->SetPoint(0, point[0], point[1], point[2]);
2722 normline->SetPoint(1, point[0] + normal[0] * normlen, point[1] + normal[1] * normlen,
2723 point[2] + normal[2] * normlen);
2724 pm->Add(normline);
2725 }
2726 ipoint = 0;
2727 line = nullptr;
2728 }
2729 if (vis2) {
2730 // create new segment
2731 line = new TPolyLine3D(2);
2732 line->SetLineColor(endnode->GetVolume()->GetLineColor());
2733 line->SetPoint(ipoint++, point[0], point[1], point[2]);
2734 i++;
2735 if (check_norm) {
2736 normline = new TPolyLine3D(2);
2737 normline->SetLineColor(kBlue);
2738 normline->SetLineWidth(2);
2739 normline->SetPoint(0, point[0], point[1], point[2]);
2740 normline->SetPoint(1, point[0] + normal[0] * normlen, point[1] + normal[1] * normlen,
2741 point[2] + normal[2] * normlen);
2742 }
2743 pm->Add(line);
2744 if (!random)
2745 pm->Add(normline);
2746 }
2747 }
2748 }
2749 // draw all segments
2750 for (Int_t m = 0; m < pm->GetEntriesFast(); m++) {
2751 line = (TPolyLine3D *)pm->At(m);
2752 if (line)
2753 line->Draw("SAME");
2754 }
2755 printf("number of segments : %i\n", i);
2756 // fGeoManager->GetTopVolume()->VisibleDaughters(kFALSE);
2757 // printf("---Daughters of %s made invisible.\n", fGeoManager->GetTopVolume()->GetName());
2758 // printf("---Make them visible with : gGeoManager->GetTopVolume()->VisibleDaughters();\n");
2759 delete pm;
2760}
2761
2762////////////////////////////////////////////////////////////////////////////////
2763/// shoot npoints randomly in a box of 1E-5 around current point.
2764/// return minimum distance to points outside
2765/// make sure that path to current node is updated
2766/// get the response of tgeo
2767
2769{
2770 TGeoNode *node = fGeoManager->FindNode();
2771 TGeoNode *nodegeo = nullptr;
2772 TGeoNode *nodeg3 = nullptr;
2773 TGeoNode *solg3 = nullptr;
2774 if (!node) {
2775 dist = -1;
2776 return nullptr;
2777 }
2779 if (strlen(g3path))
2780 hasg3 = kTRUE;
2782 dist = 1E10;
2783 TString common = "";
2784 // cd to common path
2785 Double_t point[3];
2786 Double_t closest[3];
2787 TGeoNode *node1 = nullptr;
2788 TGeoNode *node_close = nullptr;
2789 dist = 1E10;
2790 Double_t dist1 = 0;
2791 // initialize size of random box to epsil
2792 Double_t eps[3];
2793 eps[0] = epsil;
2794 eps[1] = epsil;
2795 eps[2] = epsil;
2797 if (hasg3) {
2799 TString name = "";
2800 Int_t index = 0;
2801 while (index >= 0) {
2802 index = spath.Index("/", index + 1);
2803 if (index > 0) {
2804 name = spath(0, index);
2805 if (strstr(g3path, name.Data())) {
2806 common = name;
2807 continue;
2808 } else
2809 break;
2810 }
2811 }
2812 // if g3 response was given, cd to common path
2813 if (strlen(common.Data())) {
2814 while (strcmp(fGeoManager->GetPath(), common.Data()) && fGeoManager->GetLevel()) {
2816 fGeoManager->CdUp();
2817 }
2820 while (strcmp(fGeoManager->GetPath(), common.Data()) && fGeoManager->GetLevel()) {
2822 fGeoManager->CdUp();
2823 }
2824 if (!nodegeo)
2825 return nullptr;
2826 if (!nodeg3)
2827 return nullptr;
2828 fGeoManager->cd(common.Data());
2830 Double_t xyz[3], local[3];
2831 for (Int_t i = 0; i < npoints; i++) {
2832 xyz[0] = point[0] - eps[0] + 2 * eps[0] * gRandom->Rndm();
2833 xyz[1] = point[1] - eps[1] + 2 * eps[1] * gRandom->Rndm();
2834 xyz[2] = point[2] - eps[2] + 2 * eps[2] * gRandom->Rndm();
2835 nodeg3->MasterToLocal(&xyz[0], &local[0]);
2836 if (!nodeg3->GetVolume()->Contains(&local[0]))
2837 continue;
2838 dist1 = TMath::Sqrt((xyz[0] - point[0]) * (xyz[0] - point[0]) + (xyz[1] - point[1]) * (xyz[1] - point[1]) +
2839 (xyz[2] - point[2]) * (xyz[2] - point[2]));
2840 if (dist1 < dist) {
2841 // save node and closest point
2842 dist = dist1;
2843 node_close = solg3;
2844 // make the random box smaller
2845 eps[0] = TMath::Abs(point[0] - pointg[0]);
2846 eps[1] = TMath::Abs(point[1] - pointg[1]);
2847 eps[2] = TMath::Abs(point[2] - pointg[2]);
2848 }
2849 }
2850 }
2851 if (!node_close)
2852 dist = -1;
2853 return node_close;
2854 }
2855
2856 // save current point
2857 memcpy(&point[0], pointg, 3 * sizeof(Double_t));
2858 for (Int_t i = 0; i < npoints; i++) {
2859 // generate a random point in MARS
2860 fGeoManager->SetCurrentPoint(point[0] - eps[0] + 2 * eps[0] * gRandom->Rndm(),
2861 point[1] - eps[1] + 2 * eps[1] * gRandom->Rndm(),
2862 point[2] - eps[2] + 2 * eps[2] * gRandom->Rndm());
2864 // check if new node1 is different from the old one
2865 if (node1 != node) {
2866 dist1 = TMath::Sqrt((point[0] - pointg[0]) * (point[0] - pointg[0]) +
2867 (point[1] - pointg[1]) * (point[1] - pointg[1]) +
2868 (point[2] - pointg[2]) * (point[2] - pointg[2]));
2869 if (dist1 < dist) {
2870 dist = dist1;
2871 node_close = node1;
2872 memcpy(&closest[0], pointg, 3 * sizeof(Double_t));
2873 // make the random box smaller
2874 eps[0] = TMath::Abs(point[0] - pointg[0]);
2875 eps[1] = TMath::Abs(point[1] - pointg[1]);
2876 eps[2] = TMath::Abs(point[2] - pointg[2]);
2877 }
2878 }
2879 }
2880 // restore the original point and path
2881 fGeoManager->FindNode(point[0], point[1], point[2]); // really needed ?
2882 if (!node_close)
2883 dist = -1;
2884 return node_close;
2885}
2886
2887////////////////////////////////////////////////////////////////////////////////
2888/// Shoot one ray from start point with direction (dirx,diry,dirz). Fills input array
2889/// with points just after boundary crossings.
2890
2892 Int_t &nelem, Int_t &dim, Double_t *endpoint) const
2893{
2894 // Int_t array_dimension = 3*dim;
2895 nelem = 0;
2896 Int_t istep = 0;
2897 if (!dim) {
2898 printf("empty input array\n");
2899 return array;
2900 }
2901 // fGeoManager->CdTop();
2902 const Double_t *point = fGeoManager->GetCurrentPoint();
2905 Double_t step, forward;
2906 Double_t dir[3];
2907 dir[0] = dirx;
2908 dir[1] = diry;
2909 dir[2] = dirz;
2910 fGeoManager->InitTrack(start, &dir[0]);
2912 // printf("Start : (%f,%f,%f)\n", point[0], point[1], point[2]);
2914 step = fGeoManager->GetStep();
2915 // printf("---next : at step=%f\n", step);
2916 if (step > 1E10)
2917 return array;
2920 while (step < 1E10) {
2921 if (endpoint) {
2922 forward = dirx * (endpoint[0] - point[0]) + diry * (endpoint[1] - point[1]) + dirz * (endpoint[2] - point[2]);
2923 if (forward < 1E-3) {
2924 // printf("exit : Passed start point. nelem=%i\n", nelem);
2925 return array;
2926 }
2927 }
2928 if (is_entering) {
2929 if (nelem >= dim) {
2930 Double_t *temparray = new Double_t[3 * (dim + 20)];
2931 memcpy(temparray, array, 3 * dim * sizeof(Double_t));
2932 delete[] array;
2933 array = temparray;
2934 dim += 20;
2935 }
2936 memcpy(&array[3 * nelem], point, 3 * sizeof(Double_t));
2937 // printf("%i (%f, %f, %f) step=%f\n", nelem, point[0], point[1], point[2], step);
2938 nelem++;
2939 } else {
2940 if (endnode == nullptr && step > 1E10) {
2941 // printf("exit : NULL endnode. nelem=%i\n", nelem);
2942 return array;
2943 }
2944 if (!fGeoManager->IsEntering()) {
2945 // if (startnode) printf("stepping %f from (%f, %f, %f) inside %s...\n", step,point[0], point[1],
2946 // point[2], startnode->GetName()); else printf("stepping %f from (%f, %f, %f) OUTSIDE...\n",
2947 // step,point[0], point[1], point[2]);
2948 istep = 0;
2949 }
2950 while (!fGeoManager->IsEntering()) {
2951 istep++;
2952 if (istep > 1E3) {
2953 // Error("ShootRay", "more than 1000 steps. Step was %f", step);
2954 nelem = 0;
2955 return array;
2956 }
2957 fGeoManager->SetStep(1E-5);
2959 }
2960 if (istep > 0)
2961 printf("%i steps\n", istep);
2962 if (nelem >= dim) {
2963 Double_t *temparray = new Double_t[3 * (dim + 20)];
2964 memcpy(temparray, array, 3 * dim * sizeof(Double_t));
2965 delete[] array;
2966 array = temparray;
2967 dim += 20;
2968 }
2969 memcpy(&array[3 * nelem], point, 3 * sizeof(Double_t));
2970 // if (endnode) printf("%i (%f, %f, %f) step=%f\n", nelem, point[0], point[1], point[2], step);
2971 nelem++;
2972 }
2974 step = fGeoManager->GetStep();
2975 // printf("---next at step=%f\n", step);
2978 }
2979 return array;
2980 // printf("exit : INFINITE step. nelem=%i\n", nelem);
2981}
2982
2983////////////////////////////////////////////////////////////////////////////////
2984/// Check time of finding "Where am I" for n points.
2985
2987{
2988 Bool_t recheck = !strcmp(option, "RECHECK");
2989 if (recheck)
2990 printf("RECHECK\n");
2991 const TGeoShape *shape = fGeoManager->GetTopVolume()->GetShape();
2992 Double_t dx = ((TGeoBBox *)shape)->GetDX();
2993 Double_t dy = ((TGeoBBox *)shape)->GetDY();
2994 Double_t dz = ((TGeoBBox *)shape)->GetDZ();
2995 Double_t ox = (((TGeoBBox *)shape)->GetOrigin())[0];
2996 Double_t oy = (((TGeoBBox *)shape)->GetOrigin())[1];
2997 Double_t oz = (((TGeoBBox *)shape)->GetOrigin())[2];
2998 Double_t *xyz = new Double_t[3 * npoints];
2999 TStopwatch *timer = new TStopwatch();
3000 printf("Random box : %f, %f, %f\n", dx, dy, dz);
3001 timer->Start(kFALSE);
3002 Int_t i;
3003 for (i = 0; i < npoints; i++) {
3004 xyz[3 * i] = ox - dx + 2 * dx * gRandom->Rndm();
3005 xyz[3 * i + 1] = oy - dy + 2 * dy * gRandom->Rndm();
3006 xyz[3 * i + 2] = oz - dz + 2 * dz * gRandom->Rndm();
3007 }
3008 timer->Stop();
3009 printf("Generation time :\n");
3010 timer->Print();
3011 timer->Reset();
3012 TGeoNode *node, *node1;
3013 printf("Start... %i points\n", npoints);
3014 timer->Start(kFALSE);
3015 for (i = 0; i < npoints; i++) {
3016 fGeoManager->SetCurrentPoint(xyz + 3 * i);
3017 if (recheck)
3018 fGeoManager->CdTop();
3019 node = fGeoManager->FindNode();
3020 if (recheck) {
3022 if (node1 != node) {
3023 printf("Difference for x=%g y=%g z=%g\n", xyz[3 * i], xyz[3 * i + 1], xyz[3 * i + 2]);
3024 printf(" from top : %s\n", node->GetName());
3025 printf(" redo : %s\n", fGeoManager->GetPath());
3026 }
3027 }
3028 }
3029 timer->Stop();
3030 timer->Print();
3031 delete[] xyz;
3032 delete timer;
3033}
3034
3035////////////////////////////////////////////////////////////////////////////////
3036/// Geometry overlap checker based on sampling.
3037
3038void TGeoChecker::TestOverlaps(const char *path)
3039{
3042 printf("Checking overlaps for path :\n");
3043 if (!fGeoManager->cd(path))
3044 return;
3045 TGeoNode *checked = fGeoManager->GetCurrentNode();
3046 checked->InspectNode();
3047 // shoot 1E4 points in the shape of the current volume
3048 Int_t npoints = 1000000;
3049 Double_t big = 1E6;
3050 Double_t xmin = big;
3051 Double_t xmax = -big;
3052 Double_t ymin = big;
3053 Double_t ymax = -big;
3054 Double_t zmin = big;
3055 Double_t zmax = -big;
3056 TObjArray *pm = new TObjArray(128);
3057 TPolyMarker3D *marker = nullptr;
3059 markthis->SetMarkerColor(5);
3060 TNtuple *ntpl = new TNtuple("ntpl", "random points", "x:y:z");
3062 Double_t *point = new Double_t[3];
3063 Double_t dx = ((TGeoBBox *)shape)->GetDX();
3064 Double_t dy = ((TGeoBBox *)shape)->GetDY();
3065 Double_t dz = ((TGeoBBox *)shape)->GetDZ();
3066 Double_t ox = (((TGeoBBox *)shape)->GetOrigin())[0];
3067 Double_t oy = (((TGeoBBox *)shape)->GetOrigin())[1];
3068 Double_t oz = (((TGeoBBox *)shape)->GetOrigin())[2];
3069 Double_t *xyz = new Double_t[3 * npoints];
3070 Int_t i = 0;
3071 printf("Generating %i points inside %s\n", npoints, fGeoManager->GetPath());
3072 while (i < npoints) {
3073 point[0] = ox - dx + 2 * dx * gRandom->Rndm();
3074 point[1] = oy - dy + 2 * dy * gRandom->Rndm();
3075 point[2] = oz - dz + 2 * dz * gRandom->Rndm();
3076 if (!shape->Contains(point))
3077 continue;
3078 // convert each point to MARS
3079 // printf("local %9.3f %9.3f %9.3f\n", point[0], point[1], point[2]);
3080 fGeoManager->LocalToMaster(point, &xyz[3 * i]);
3081 // printf("master %9.3f %9.3f %9.3f\n", xyz[3*i], xyz[3*i+1], xyz[3*i+2]);
3082 xmin = TMath::Min(xmin, xyz[3 * i]);
3083 xmax = TMath::Max(xmax, xyz[3 * i]);
3084 ymin = TMath::Min(ymin, xyz[3 * i + 1]);
3085 ymax = TMath::Max(ymax, xyz[3 * i + 1]);
3086 zmin = TMath::Min(zmin, xyz[3 * i + 2]);
3087 zmax = TMath::Max(zmax, xyz[3 * i + 2]);
3088 i++;
3089 }
3090 delete[] point;
3091 ntpl->Fill(xmin, ymin, zmin);
3092 ntpl->Fill(xmax, ymin, zmin);
3093 ntpl->Fill(xmin, ymax, zmin);
3094 ntpl->Fill(xmax, ymax, zmin);
3095 ntpl->Fill(xmin, ymin, zmax);
3096 ntpl->Fill(xmax, ymin, zmax);
3097 ntpl->Fill(xmin, ymax, zmax);
3098 ntpl->Fill(xmax, ymax, zmax);
3099 ntpl->Draw("z:y:x");
3100
3101 // shoot the poins in the geometry
3102 TGeoNode *node;
3103 TString cpath;
3104 Int_t ic = 0;
3105 TObjArray *overlaps = new TObjArray();
3106 printf("using FindNode...\n");
3107 for (Int_t j = 0; j < npoints; j++) {
3108 // always start from top level (testing only)
3109 fGeoManager->CdTop();
3110 fGeoManager->SetCurrentPoint(&xyz[3 * j]);
3111 node = fGeoManager->FindNode();
3113 if (cpath.Contains(path)) {
3114 markthis->SetNextPoint(xyz[3 * j], xyz[3 * j + 1], xyz[3 * j + 2]);
3115 continue;
3116 }
3117 // current point is found in an overlapping node
3118 if (!node)
3119 ic = 128;
3120 else
3121 ic = node->GetVolume()->GetLineColor();
3122 if (ic >= 128)
3123 ic = 0;
3124 marker = (TPolyMarker3D *)pm->At(ic);
3125 if (!marker) {
3126 marker = new TPolyMarker3D();
3127 marker->SetMarkerColor(ic);
3128 pm->AddAt(marker, ic);
3129 }
3130 // draw the overlapping point
3131 marker->SetNextPoint(xyz[3 * j], xyz[3 * j + 1], xyz[3 * j + 2]);
3132 if (node) {
3133 if (overlaps->IndexOf(node) < 0)
3134 overlaps->Add(node);
3135 }
3136 }
3137 // draw all overlapping points
3138 // for (Int_t m=0; m<128; m++) {
3139 // marker = (TPolyMarker3D*)pm->At(m);
3140 // if (marker) marker->Draw("SAME");
3141 // }
3142 markthis->Draw("SAME");
3143 if (gPad)
3144 gPad->Update();
3145 // display overlaps
3146 if (overlaps->GetEntriesFast()) {
3147 printf("list of overlapping nodes :\n");
3148 for (i = 0; i < overlaps->GetEntriesFast(); i++) {
3149 node = (TGeoNode *)overlaps->At(i);
3150 if (node->IsOverlapping())
3151 printf("%s MANY\n", node->GetName());
3152 else
3153 printf("%s ONLY\n", node->GetName());
3154 }
3155 } else
3156 printf("No overlaps\n");
3157 delete ntpl;
3158 delete pm;
3159 delete[] xyz;
3160 delete overlaps;
3161}
3162
3163////////////////////////////////////////////////////////////////////////////////
3164/// Estimate weight of top level volume with a precision SIGMA(W)/W
3165/// better than PRECISION. Option can be "v" - verbose (default).
3166
3168{
3170 Int_t nmat = matlist->GetSize();
3171 if (!nmat)
3172 return 0;
3173 Int_t *nin = new Int_t[nmat];
3174 memset(nin, 0, nmat * sizeof(Int_t));
3175 TString opt = option;
3176 opt.ToLower();
3177 Bool_t isverbose = opt.Contains("v");
3179 Double_t dx = box->GetDX();
3180 Double_t dy = box->GetDY();
3181 Double_t dz = box->GetDZ();
3182 Double_t ox = (box->GetOrigin())[0];
3183 Double_t oy = (box->GetOrigin())[1];
3184 Double_t oz = (box->GetOrigin())[2];
3185 Double_t x, y, z;
3186 TGeoNode *node;
3188 Double_t vbox = 0.000008 * dx * dy * dz; // m3
3189 Bool_t end = kFALSE;
3190 Double_t weight = 0, sigma, eps, dens;
3191 Double_t eps0 = 1.;
3192 Int_t indmat;
3193 Int_t igen = 0;
3194 Int_t iin = 0;
3195 while (!end) {
3196 x = ox - dx + 2 * dx * gRandom->Rndm();
3197 y = oy - dy + 2 * dy * gRandom->Rndm();
3198 z = oz - dz + 2 * dz * gRandom->Rndm();
3199 node = fGeoManager->FindNode(x, y, z);
3200 igen++;
3201 if (!node)
3202 continue;
3203 mat = node->GetVolume()->GetMedium()->GetMaterial();
3204 indmat = mat->GetIndex();
3205 if (indmat < 0)
3206 continue;
3207 nin[indmat]++;
3208 iin++;
3209 if ((iin % 100000) == 0 || igen > 1E8) {
3210 weight = 0;
3211 sigma = 0;
3212 for (indmat = 0; indmat < nmat; indmat++) {
3213 mat = (TGeoMaterial *)matlist->At(indmat);
3214 dens = mat->GetDensity(); // [g/cm3]
3215 if (dens < 1E-2)
3216 dens = 0;
3217 dens *= 1000.; // [kg/m3]
3218 weight += dens * Double_t(nin[indmat]);
3219 sigma += dens * dens * nin[indmat];
3220 }
3222 eps = sigma / weight;
3223 weight *= vbox / Double_t(igen);
3224 sigma *= vbox / Double_t(igen);
3226 if (isverbose) {
3227 printf("=== Weight of %s : %g +/- %g [kg]\n", fGeoManager->GetTopVolume()->GetName(), weight, sigma);
3228 }
3229 end = kTRUE;
3230 } else {
3231 if (isverbose && eps < 0.5 * eps0) {
3232 printf("%8dK: %14.7g kg %g %%\n", igen / 1000, weight, eps * 100);
3233 eps0 = eps;
3234 }
3235 }
3236 }
3237 }
3238 delete[] nin;
3239 return weight;
3240}
3241////////////////////////////////////////////////////////////////////////////////
3242/// count voxel timing
3243
3245{
3247 Double_t time;
3248 TGeoShape *shape = vol->GetShape();
3249 TGeoNode *node;
3251 Double_t *point;
3252 Double_t local[3];
3254 Int_t ncheck;
3256 TGeoStateInfo &td = *nav->GetCache()->GetInfo();
3257 timer.Start();
3258 for (Int_t i = 0; i < npoints; i++) {
3259 point = xyz + 3 * i;
3260 if (!shape->Contains(point))
3261 continue;
3262 checklist = voxels->GetCheckList(point, ncheck, td);
3263 if (!checklist)
3264 continue;
3265 if (!ncheck)
3266 continue;
3267 for (Int_t id = 0; id < ncheck; id++) {
3268 node = vol->GetNode(checklist[id]);
3269 matrix = node->GetMatrix();
3270 matrix->MasterToLocal(point, &local[0]);
3271 if (node->GetVolume()->GetShape()->Contains(&local[0]))
3272 break;
3273 }
3274 }
3275 nav->GetCache()->ReleaseInfo();
3276 time = timer.CpuTime();
3277 return time;
3278}
3279
3280////////////////////////////////////////////////////////////////////////////////
3281/// Returns optimal voxelization type for volume vol.
3282/// - kFALSE - cartesian
3283/// - kTRUE - cylindrical
3284
3286{
3287 return kFALSE;
3288}
std::ios_base::fmtflags 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
bool Bool_t
Boolean (0=false, 1=true) (bool)
Definition RtypesCore.h:77
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
char Char_t
Character 1 byte (char)
Definition RtypesCore.h:51
long Long_t
Signed long integer 4 bytes (long). Size depends on architecture.
Definition RtypesCore.h:68
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
long long Long64_t
Portable signed long integer 8 bytes.
Definition RtypesCore.h:83
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
@ kRed
Definition Rtypes.h:67
@ kBlue
Definition Rtypes.h:67
@ kYellow
Definition Rtypes.h:67
@ kFullCircle
Definition TAttMarker.h:58
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
winID h TVirtualViewer3D TVirtualGLPainter p
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 TGeoManager * gGeoManager
R__EXTERN TGeoIdentity * gGeoIdentity
Definition TGeoMatrix.h:538
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:2495
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
#define gPad
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:38
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:35
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:42
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:39
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:41
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition TAxis.cxx:482
Generic 3D primitive description class.
Definition TBuffer3D.h:18
The Canvas class.
Definition TCanvas.h:23
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:130
Box class.
Definition TGeoBBox.h:18
static bool MayIntersect(const TGeoBBox *boxA, const TGeoMatrix *mA, const TGeoBBox *boxB, const TGeoMatrix *mB)
Fast "may-intersect" test for two placed TGeoBBox objects.
Definition TGeoBBox.cxx:952
Int_t FillMeshPoints(TBuffer3D &buff, const TGeoShape *shape, const TGeoShape *&lastShape, Int_t &cachedN, Int_t nMeshPoints, Double_t *&points) const
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:66
Bool_t * fFlags
Array of timing per volume.
Definition TGeoChecker.h:63
TStopwatch * fTimer
Array of flags per volume.
Definition TGeoChecker.h:64
void CheckShape(TGeoShape *shape, Int_t testNo, Int_t nsamples, Option_t *option) override
Test for shape navigation methods.
TGeoVolume * fVsafe
Definition TGeoChecker.h:59
std::unordered_map< const TGeoShape *, MeshPointsEntry > fMeshPointsCache
Definition TGeoChecker.h:69
TGeoChecker()
Default constructor.
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="") override
Print current operation progress.
TGeoNode * SamplePoints(Int_t npoints, Double_t &dist, Double_t epsil, const char *g3path) override
shoot npoints randomly in a box of 1E-5 around current point.
void CheckPoint(Double_t x=0, Double_t y=0, Double_t z=0, Option_t *option="", Double_t safety=0.) override
Draw point (x,y,z) over the picture of the daughters of the volume containing this point.
Double_t Weight(Double_t precision=0.01, Option_t *option="v") override
Estimate weight of top level volume with a precision SIGMA(W)/W better than PRECISION.
void ShapeNormal(TGeoShape *shape, Int_t nsamples, Option_t *option)
Check of validity of the normal for a given shape.
MeshPolicyStamp fMeshPolicyStamp
Internal map shape->mesh points.
Definition TGeoChecker.h:70
void CheckBoundaryReference(Int_t icheck=-1) override
Check the boundary errors reference file created by CheckBoundaryErrors method.
void MaterializeOverlap(const TGeoOverlapResult &r) override
Materialize a TGeoOverlapResult into a TGeoOverlap.
Double_t * fVal2
Array of number of crossings per volume.
Definition TGeoChecker.h:62
const MeshPointsEntry & GetMeshPoints(const TGeoShape *s) const
Bool_t TestVoxels(TGeoVolume *vol, Int_t npoints=1000000) override
Returns optimal voxelization type for volume vol.
void CheckOverlapsBySampling(const TGeoVolume *vol, Double_t ovlp=0.1, Int_t npoints=1000000) const override
Check illegal overlaps for volume VOL within a limit OVLP by sampling npoints inside the volume shape...
void BuildMeshPointsCache(const std::vector< TGeoOverlapCandidate > &candidates) override
TGeoManager * fGeoManager
Definition TGeoChecker.h:58
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.
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="") override
Generate a lego plot fot the top volume, according to option.
void CheckGeometry(Int_t nrays, Double_t startx, Double_t starty, Double_t startz) const override
Shoot nrays with random directions from starting point (startx, starty, startz) in the reference fram...
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
policy stamp for points cache invalidation
bool PolicyChanged(const MeshPolicyStamp &a, const MeshPolicyStamp &b) const
Definition TGeoChecker.h:79
void PrintOverlaps() const override
Print the current list of overlaps held by the manager class.
TGeoNode * fSelectedNode
Timer.
Definition TGeoChecker.h:65
void CheckGeometryFull(Bool_t checkoverlaps=kTRUE, Bool_t checkcrossings=kTRUE, Int_t nrays=10000, const Double_t *vertex=nullptr) override
Geometry checking.
void CheckBoundaryErrors(Int_t ntracks=1000000, Double_t radius=-1.) override
Check pushes and pulls needed to cross the next boundary with respect to the position given by FindNe...
~TGeoChecker() override
Destructor.
Bool_t ComputeOverlap(const TGeoOverlapCandidate &c, TGeoOverlapResult &out) const override
Compute overlap/extrusion for given candidate using mesh points of the shapes.
TGeoChecker::MeshPolicyStamp CurrentMeshPolicyStamp() const
Create a policy stamp for overlap checking.
void TestOverlaps(const char *path) override
Geometry overlap checker based on sampling.
void Test(Int_t npoints, Option_t *option) override
Check time of finding "Where am I" for n points.
std::shared_mutex fMeshCacheMutex
Number of checks for current volume.
Definition TGeoChecker.h:68
Int_t EnumerateOverlapCandidates(const TGeoVolume *vol, Double_t ovlp, Option_t *option, std::vector< TGeoOverlapCandidate > &out) override
Stage1: Enumerate all overlap candidates for volume VOL within a limit OVLP.
Double_t TimingPerVolume(TGeoVolume *)
Compute timing per "FindNextBoundary" + "Safety" call.
void RandomRays(Int_t nrays, Double_t startx, Double_t starty, Double_t startz, const char *target_vol=nullptr, Bool_t check_norm=kFALSE) override
Randomly shoot nrays from point (startx,starty,startz) and plot intersections with surfaces for curre...
Double_t * fVal1
Definition TGeoChecker.h:61
Bool_t fFullCheck
Definition TGeoChecker.h:60
void RandomPoints(TGeoVolume *vol, Int_t npoints, Option_t *option) override
Draw random points in the bounding box of a volume.
Double_t CheckVoxels(TGeoVolume *vol, TGeoVoxelFinder *voxels, Double_t *xyz, Int_t npoints)
count voxel timing
void PushCandidate(std::vector< TGeoOverlapCandidate > &out, const TString &name, TGeoVolume *vol1, TGeoVolume *vol2, const TGeoMatrix *mat1, const TGeoMatrix *mat2, Bool_t isovlp, Double_t ovlp) const
Helper to fill candidates list.
Matrix class used for computing global transformations Should NOT be used for node definition.
Definition TGeoMatrix.h:459
A geometry iterator.
Definition TGeoNode.h:249
void GetPath(TString &path) const
Returns the path for the current node.
The manager class for any TGeo geometry.
Definition TGeoManager.h:46
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
Int_t GetNsegments() const
Get number of segments approximating circles.
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.
Geometrical transformation package.
Definition TGeoMatrix.h:39
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 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:49
Class providing navigation API for TGeo geometries.
A node represent a volume positioned inside another.They store links to both volumes and to the TGeoM...
Definition TGeoNode.h:39
Bool_t IsOverlapping() const
Definition TGeoNode.h:108
Bool_t IsOnScreen() const
check if this node is drawn. Assumes that this node is current
Definition TGeoNode.cxx:422
TGeoVolume * GetVolume() const
Definition TGeoNode.h:100
virtual TGeoMatrix * GetMatrix() const =0
Int_t GetColour() const
Definition TGeoNode.h:87
void InspectNode() const
Inspect this node.
Definition TGeoNode.cxx:432
Base class describing geometry overlaps.
Definition TGeoOverlap.h:37
Base abstract class for all shapes.
Definition TGeoShape.h:25
static Double_t Big()
Definition TGeoShape.h:95
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:133
virtual Bool_t GetPointsOnSegments(Int_t npoints, Double_t *array) const =0
void Draw(Option_t *option="") override
Draw this shape.
virtual Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const =0
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm) const =0
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
const char * GetName() const override
Get the shape name.
virtual Double_t Capacity() const =0
virtual Bool_t Contains(const Double_t *point) const =0
static Double_t Tolerance()
Definition TGeoShape.h:98
virtual void SetPoints(Double_t *points) const =0
Tessellated solid class.
Class describing translations.
Definition TGeoMatrix.h:117
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:43
void SetVisibility(Bool_t vis=kTRUE) override
set visibility of this volume
TGeoMedium * GetMedium() const
Definition TGeoVolume.h:176
Bool_t Contains(const Double_t *point) const
Definition TGeoVolume.h:105
virtual TGeoNode * AddNode(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat=nullptr, Option_t *option="")
Add a TGeoNode to the list of nodes.
void Draw(Option_t *option="") override
draw top volume according to option
Int_t GetNdaughters() const
Definition TGeoVolume.h:363
TObjArray * GetNodes()
Definition TGeoVolume.h:170
void VisibleDaughters(Bool_t vis=kTRUE)
set visibility for daughters
TGeoNode * GetNode(const char *name) const
get the pointer to a daughter node
Int_t GetIndex(const TGeoNode *node) const
get index number for a given daughter
TGeoPatternFinder * GetFinder() const
Definition TGeoVolume.h:178
TGeoVoxelFinder * GetVoxels() const
Getter for optimization structure.
Int_t GetNumber() const
Definition TGeoVolume.h:185
TGeoShape * GetShape() const
Definition TGeoVolume.h:191
virtual void DrawOnly(Option_t *option="")
draw only this volume
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:156
void InspectShape() const
Definition TGeoVolume.h:196
Finder class handling voxels.
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:926
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:878
@ kAllAxes
Definition TH1.h:126
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:5402
static void AddDirectory(Bool_t add=kTRUE)
Sets the flag controlling the automatic add of histograms in memory.
Definition TH1.cxx:1263
TAxis * GetXaxis()
Definition TH1.h:571
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3330
TAxis * GetYaxis()
Definition TH1.h:572
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3052
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:6702
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:5265
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition TH1.cxx:9050
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:345
Int_t Fill(Double_t) override
Invalid Fill method.
Definition TH2.cxx:363
A doubly linked list.
Definition TList.h:38
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:173
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
A simple TTree restricted to a list of float variables only.
Definition TNtuple.h:28
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
Int_t IndexOf(const TObject *obj) const override
Int_t GetEntries() const override
Return the number of objects in array (i.e.
TObject * At(Int_t idx) const override
Definition TObjArray.h:170
TObject * RemoveAt(Int_t idx) override
Remove object at index idx.
void Add(TObject *obj) override
Definition TObjArray.h:68
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:226
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1074
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1088
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:293
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:1062
A 3-dimensional polyline.
Definition TPolyLine3D.h:33
A 3D polymarker.
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:558
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition TRandom.cxx:681
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 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:138
void ToLower()
Change string to lower-case.
Definition TString.cxx:1189
const char * Data() const
Definition TString.h:384
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:2384
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:641
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:1641
A TTree represents a columnar dataset.
Definition TTree.h:89
Abstract class for geometry checkers.
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:643
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:249
Double_t Log(Double_t x)
Returns the natural logarithm of x.
Definition TMath.h:767
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:673
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:197
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:605
constexpr Double_t Pi()
Definition TMath.h:40
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:599
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:122
constexpr Double_t TwoPi()
Definition TMath.h:47
Simple internal structure keeping a cache of mesh points per shape.
Definition TGeoChecker.h:45
Policy stamp: when this changes, cached point values may be wrong -> clear.
Definition TGeoChecker.h:51
TString fName
display name
Statefull info for the current geometry level.
TMarker m
Definition textangle.C:8