Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoVolume.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Andrei Gheata 30/05/02
3// Divide(), CheckOverlaps() implemented 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 TGeoVolume
14\ingroup Shapes_classes
15
16TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes
17
18 Volumes are the basic objects used in building the geometrical hierarchy.
19They represent unpositioned objects but store all information about the
20placement of the other volumes they may contain. Therefore a volume can
21be replicated several times in the geometry. In order to create a volume, one
22has to put together a shape and a medium which are already defined. Volumes
23have to be named by users at creation time. Every different name may represent a
24an unique volume object, but may also represent more general a family (class)
25of volume objects having the same shape type and medium, but possibly
26different shape parameters. It is the user's task to provide different names
27for different volume families in order to avoid ambiguities at tracking time.
28A generic family rather than a single volume is created only in two cases :
29when a generic shape is provided to the volume constructor or when a division
30operation is applied. Each volume in the geometry stores an unique
31ID corresponding to its family. In order to ease-up their creation, the manager
32class is providing an API that allows making a shape and a volume in a single step.
33
34 Volumes are objects that can be visualized, therefore having visibility,
35colour, line and fill attributes that can be defined or modified any time after
36the volume creation. It is advisable however to define these properties just
37after the first creation of a volume namespace, since in case of volume families
38any new member created by the modeler inherits these properties.
39
40 In order to provide navigation features, volumes have to be able to find
41the proper container of any point defined in the local reference frame. This
42can be the volume itself, one of its positioned daughter volumes or none if
43the point is actually outside. On the other hand, volumes have to provide also
44other navigation methods such as finding the distances to its shape boundaries
45or which daughter will be crossed first. The implementation of these features
46is done at shape level, but the local mother-daughters management is handled
47by volumes that builds additional optimisation structures upon geometry closure.
48In order to have navigation features properly working one has to follow the
49general rules for building a valid geometry (see TGeoManager class).
50
51 Now let's make a simple volume representing a copper wire. We suppose that
52a medium is already created (see TGeoMedium class on how to create media).
53We will create a TUBE shape for our wire, having Rmin=0cm, Rmax=0.01cm
54and a half-length dZ=1cm :
55
56~~~ {.cpp}
57 TGeoTube *tube = new TGeoTube("wire_tube", 0, 0.01, 1);
58~~~
59
60One may omit the name for the shape if no retrieving by name is further needed
61during geometry building. The same shape can be shared by different volumes
62having different names and materials. Now let's make the volume for our wire.
63The prototype for volumes constructor looks like :
64
65 TGeoVolume::TGeoVolume(const char *name, TGeoShape *shape, TGeoMedium *med)
66
67Since TGeoTube derives from the base shape class, we can provide it to the volume
68constructor :
69
70~~~ {.cpp}
71 TGeoVolume *wire_co = new TGeoVolume("WIRE_CO", tube, ptrCOPPER);
72~~~
73
74Do not bother to delete neither the media, shapes or volumes that you have
75created since all will be automatically cleaned on exit by the manager class.
76If we would have taken a look inside TGeoManager::MakeTube() method, we would
77have been able to create our wire with a single line :
78
79~~~ {.cpp}
80 TGeoVolume *wire_co = gGeoManager->MakeTube("WIRE_CO", ptrCOPPER, 0, 0.01, 1);
81~~~
82
83The same applies for all primitive shapes, for which there can be found
84corresponding MakeSHAPE() methods. Their usage is much more convenient unless
85a shape has to be shared between more volumes. Let's make now an aluminium wire
86having the same shape, supposing that we have created the copper wire with the
87line above :
88
89~~~ {.cpp}
90 TGeoVolume *wire_al = new TGeoVolume("WIRE_AL", wire_co->GetShape(), ptrAL);
91~~~
92
93Now that we have learned how to create elementary volumes, let's see how we
94can create a geometrical hierarchy.
95
96
97### Positioning volumes
98
99 When creating a volume one does not specify if this will contain or not other
100volumes. Adding daughters to a volume implies creating those and adding them
101one by one to the list of daughters. Since the volume has to know the position
102of all its daughters, we will have to supply at the same time a geometrical
103transformation with respect to its local reference frame for each of them.
104The objects referencing a volume and a transformation are called NODES and
105their creation is fully handled by the modeler. They represent the link
106elements in the hierarchy of volumes. Nodes are unique and distinct geometrical
107objects ONLY from their container point of view. Since volumes can be replicated
108in the geometry, the same node may be found on different branches.
109
110\image html geom_t_example.png width=600px
111
112 An important observation is that volume objects are owned by the TGeoManager
113class. This stores a list of all volumes in the geometry, that is cleaned
114upon destruction.
115
116 Let's consider positioning now our wire in the middle of a gas chamber. We
117need first to define the gas chamber :
118
119~~~ {.cpp}
120 TGeoVolume *chamber = gGeoManager->MakeTube("CHAMBER", ptrGAS, 0, 1, 1);
121~~~
122
123Now we can put the wire inside :
124
125~~~ {.cpp}
126 chamber->AddNode(wire_co, 1);
127~~~
128
129If we inspect now the chamber volume in a browser, we will notice that it has
130one daughter. Of course the gas has some container also, but let's keep it like
131that for the sake of simplicity. The full prototype of AddNode() is :
132
133~~~ {.cpp}
134 TGeoVolume::AddNode(TGeoVolume *daughter, Int_t usernumber,
135 TGeoMatrix *matrix=gGeoIdentity)
136~~~
137
138Since we did not supplied the third argument, the wire will be positioned with
139an identity transformation inside the chamber. One will notice that the inner
140radii of the wire and chamber are both zero - therefore, aren't the two volumes
141overlapping ? The answer is no, the modeler is even relaying on the fact that
142any daughter is fully contained by its mother. On the other hand, neither of
143the nodes positioned inside a volume should overlap with each other. We will
144see that there are allowed some exceptions to those rules.
145
146### Overlapping volumes
147
148 Positioning volumes that does not overlap their neighbours nor extrude
149their container is sometimes quite strong constraint. Some parts of the geometry
150might overlap naturally, e.g. two crossing tubes. The modeller supports such
151cases only if the overlapping nodes are declared by the user. In order to do
152that, one should use TGeoVolume::AddNodeOverlap() instead of TGeoVolume::AddNode().
153 When 2 or more positioned volumes are overlapping, not all of them have to
154be declared so, but at least one. A point inside an overlapping region equally
155belongs to all overlapping nodes, but the way these are defined can enforce
156the modeler to give priorities.
157 The general rule is that the deepest node in the hierarchy containing a point
158have the highest priority. For the same geometry level, non-overlapping is
159prioritised over overlapping. In order to illustrate this, we will consider
160few examples. We will designate non-overlapping nodes as ONLY and the others
161MANY as in GEANT3, where this concept was introduced:
162 1. The part of a MANY node B extruding its container A will never be "seen"
163during navigation, as if B was in fact the result of the intersection of A and B.
164 2. If we have two nodes A (ONLY) and B (MANY) inside the same container, all
165points in the overlapping region of A and B will be designated as belonging to A.
166 3. If A an B in the above case were both MANY, points in the overlapping
167part will be designated to the one defined first. Both nodes must have the
168same medium.
169 4. The slices of a divided MANY will be as well MANY.
170
171One needs to know that navigation inside geometry parts MANY nodes is much
172slower. Any overlapping part can be defined based on composite shapes - this
173is always recommended.
174
175### Replicating volumes
176
177 What can we do if our chamber contains two identical wires instead of one ?
178What if then we would need 1000 chambers in our detector ? Should we create
1792000 wires and 1000 chamber volumes ? No, we will just need to replicate the
180ones that we have already created.
181
182~~~ {.cpp}
183 chamber->AddNode(wire_co, 1, new TGeoTranslation(-0.2,0,0));
184 chamber->AddNode(wire_co, 2, new TGeoTranslation(0.2,0,0));
185~~~
186
187 The 2 nodes that we have created inside chamber will both point to a wire_co
188object, but will be completely distinct : WIRE_CO_1 and WIRE_CO_2. We will
189want now to place symmetrically 1000 chambers on a pad, following a pattern
190of 20 rows and 50 columns. One way to do this will be to replicate our chamber
191by positioning it 1000 times in different positions of the pad. Unfortunately,
192this is far from being the optimal way of doing what we want.
193Imagine that we would like to find out which of the 1000 chambers is containing
194a (x,y,z) point defined in the pad reference. You will never have to do that,
195since the modeller will take care of it for you, but let's guess what it has
196to do. The most simple algorithm will just loop over all daughters, convert
197the point from mother to local reference and check if the current chamber
198contains the point or not. This might be efficient for pads with few chambers,
199but definitely not for 1000. Fortunately the modeler is smarter than that and
200create for each volume some optimization structures called voxels (see Voxelization)
201to minimize the penalty having too many daughters, but if you have 100 pads like
202this in your geometry you will anyway loose a lot in your tracking performance.
203
204 The way out when volumes can be arranged according to simple patterns is the
205usage of divisions. We will describe them in detail later on. Let's think now
206at a different situation : instead of 1000 chambers of the same type, we may
207have several types of chambers. Let's say all chambers are cylindrical and have
208a wire inside, but their dimensions are different. However, we would like all
209to be represented by a single volume family, since they have the same properties.
210*/
211
212/** \class TGeoVolumeMulti
213\ingroup Geometry_classes
214
215Volume families
216
217A volume family is represented by the class TGeoVolumeMulti. It represents
218a class of volumes having the same shape type and each member will be
219identified by the same name and volume ID. Any operation applied to a
220TGeoVolume equally affects all volumes in that family. The creation of a
221family is generally not a user task, but can be forced in particular cases:
222
223~~~ {.cpp}
224 TGeoManager::Volume(const char *vname, const char *shape, Int_t nmed);
225~~~
226
227where VNAME is the family name, NMED is the medium number and SHAPE is the
228shape type that can be:
229
230~~~ {.cpp}
231 box - for TGeoBBox
232 trd1 - for TGeoTrd1
233 trd2 - for TGeoTrd2
234 trap - for TGeoTrap
235 gtra - for TGeoGtra
236 para - for TGeoPara
237 tube, tubs - for TGeoTube, TGeoTubeSeg
238 cone, cons - for TGeoCone, TgeoCons
239 eltu - for TGeoEltu
240 ctub - for TGeoCtub
241 pcon - for TGeoPcon
242 pgon - for TGeoPgon
243~~~
244
245Volumes are then added to a given family upon adding the generic name as node
246inside other volume:
247
248~~~ {.cpp}
249 TGeoVolume *box_family = gGeoManager->Volume("BOXES", "box", nmed);
250 ...
251 gGeoManager->Node("BOXES", Int_t copy_no, "mother_name",
252 Double_t x, Double_t y, Double_t z, Int_t rot_index,
253 Bool_t is_only, Double_t *upar, Int_t npar);
254~~~
255
256here:
257
258~~~ {.cpp}
259 BOXES - name of the family of boxes
260 copy_no - user node number for the created node
261 mother_name - name of the volume to which we want to add the node
262 x,y,z - translation components
263 rot_index - indx of a rotation matrix in the list of matrices
264 upar - array of actual shape parameters
265 npar - number of parameters
266~~~
267
268The parameters order and number are the same as in the corresponding shape
269constructors.
270
271 Another particular case where volume families are used is when we want
272that a volume positioned inside a container to match one ore more container
273limits. Suppose we want to position the same box inside 2 different volumes
274and we want the Z size to match the one of each container:
275
276~~~ {.cpp}
277 TGeoVolume *container1 = gGeoManager->MakeBox("C1", imed, 10,10,30);
278 TGeoVolume *container2 = gGeoManager->MakeBox("C2", imed, 10,10,20);
279 TGeoVolume *pvol = gGeoManager->MakeBox("PVOL", jmed, 3,3,-1);
280 container1->AddNode(pvol, 1);
281 container2->AddNode(pvol, 1);
282~~~
283
284 Note that the third parameter of PVOL is negative, which does not make sense
285as half-length on Z. This is interpreted as: when positioned, create a box
286replacing all invalid parameters with the corresponding dimensions of the
287container. This is also internally handled by the TGeoVolumeMulti class, which
288does not need to be instantiated by users.
289
290### Dividing volumes
291
292 Volumes can be divided according a pattern. The most simple division can
293be done along one axis, that can be: X, Y, Z, Phi, Rxy or Rxyz. Let's take
294the most simple case: we would like to divide a box in N equal slices along X
295coordinate, representing a new volume family. Supposing we already have created
296the initial box, this can be done like:
297
298~~~ {.cpp}
299 TGeoVolume *slicex = box->Divide("SLICEX", 1, N);
300~~~
301
302where SLICE is the name of the new family representing all slices and 1 is the
303slicing axis. The meaning of the axis index is the following: for all volumes
304having shapes like box, trd1, trd2, trap, gtra or para - 1,2,3 means X,Y,Z; for
305tube, tubs, cone, cons - 1 means Rxy, 2 means phi and 3 means Z; for pcon and
306pgon - 2 means phi and 3 means Z; for spheres 1 means R and 2 means phi.
307 In fact, the division operation has the same effect as positioning volumes
308in a given order inside the divided container - the advantage being that the
309navigation in such a structure is much faster. When a volume is divided, a
310volume family corresponding to the slices is created. In case all slices can
311be represented by a single shape, only one volume is added to the family and
312positioned N times inside the divided volume, otherwise, each slice will be
313represented by a distinct volume in the family.
314 Divisions can be also performed in a given range of one axis. For that, one
315have to specify also the starting coordinate value and the step:
316
317~~~ {.cpp}
318 TGeoVolume *slicex = box->Divide("SLICEX", 1, N, start, step);
319~~~
320
321A check is always done on the resulting division range : if not fitting into
322the container limits, an error message is posted. If we will browse the divided
323volume we will notice that it will contain N nodes starting with index 1 upto
324N. The first one has the lower X limit at START position, while the last one
325will have the upper X limit at START+N*STEP. The resulting slices cannot
326be positioned inside an other volume (they are by default positioned inside the
327divided one) but can be further divided and may contain other volumes:
328
329~~~ {.cpp}
330 TGeoVolume *slicey = slicex->Divide("SLICEY", 2, N1);
331 slicey->AddNode(other_vol, index, some_matrix);
332~~~
333
334 When doing that, we have to remember that SLICEY represents a family, therefore
335all members of the family will be divided on Y and the other volume will be
336added as node inside all.
337 In the example above all the resulting slices had the same shape as the
338divided volume (box). This is not always the case. For instance, dividing a
339volume with TUBE shape on PHI axis will create equal slices having TUBESEG
340shape. Other divisions can also create slices having shapes with different
341dimensions, e.g. the division of a TRD1 volume on Z.
342 When positioning volumes inside slices, one can do it using the generic
343volume family (e.g. slicey). This should be done as if the coordinate system
344of the generic slice was the same as the one of the divided volume. The generic
345slice in case of PHI division is centered with respect to X axis. If the
346family contains slices of different sizes, any volume positioned inside should
347fit into the smallest one.
348 Examples for specific divisions according to shape types can be found inside
349shape classes.
350
351~~~ {.cpp}
352 TGeoVolume::Divide(N, Xmin, Xmax, "X");
353~~~
354
355 The GEANT3 option MANY is supported by TGeoVolumeOverlap class. An overlapping
356volume is in fact a virtual container that does not represent a physical object.
357It contains a list of nodes that are not its daughters but that must be checked
358always before the container itself. This list must be defined by users and it
359is checked and resolved in a priority order. Note that the feature is non-standard
360to geometrical modelers and it was introduced just to support conversions of
361GEANT3 geometries, therefore its extensive usage should be avoided.
362*/
363
364/** \class TGeoVolumeAssembly
365\ingroup Geometry_classes
366
367Volume assemblies
368
369Assemblies a volumes that have neither a shape or a material/medium. Assemblies
370behave exactly like normal volumes grouping several daughters together, but
371the daughters can never extrude the assembly since this has no shape. However,
372a bounding box and a voxelization structure are built for assemblies as for
373normal volumes, so that navigation is still optimized. Assemblies are useful
374for grouping hierarchically volumes which are otherwise defined in a flat
375manner, but also to avoid clashes between container shapes.
376To define an assembly one should just input a name, then start adding other
377volumes (or volume assemblies) as content.
378*/
379
380#include <fstream>
381#include <iomanip>
382
383#include "TString.h"
384#include "TBuffer.h"
385#include "TBrowser.h"
386#include "TStyle.h"
387#include "TH2F.h"
388#include "TROOT.h"
389#include "TEnv.h"
390#include "TMap.h"
391#include "TFile.h"
392#include "TKey.h"
393
394#include "TGeoManager.h"
395#include "TGeoNode.h"
396#include "TGeoMatrix.h"
397#include "TVirtualGeoPainter.h"
398#include "TGeoVolume.h"
399#include "TGeoShapeAssembly.h"
400#include "TGeoScaledShape.h"
401#include "TGeoCompositeShape.h"
402#include "TGeoVoxelFinder.h"
403#include "TGeoExtension.h"
404
406
408
409////////////////////////////////////////////////////////////////////////////////
410/// Create a dummy medium
411
413{
414 if (fgDummyMedium) return;
416 fgDummyMedium->SetName("dummy");
417 TGeoMaterial *dummyMaterial = new TGeoMaterial();
418 dummyMaterial->SetName("dummy");
419 fgDummyMedium->SetMaterial(dummyMaterial);
420}
421
422////////////////////////////////////////////////////////////////////////////////
423
425{
428}
429
430////////////////////////////////////////////////////////////////////////////////
431
433{
434 if (fFinder) fFinder->CreateThreadData(nthreads);
435 if (fShape) fShape->CreateThreadData(nthreads);
436}
437
438////////////////////////////////////////////////////////////////////////////////
439
441{
442 return fgDummyMedium;
443}
444
445////////////////////////////////////////////////////////////////////////////////
446/// dummy constructor
447
449{
450 fNodes = 0;
451 fShape = 0;
452 fMedium = 0;
453 fFinder = 0;
454 fVoxels = 0;
456 fField = 0;
457 fOption = "";
458 fNumber = 0;
459 fNtotal = 0;
460 fRefCount = 0;
461 fUserExtension = 0;
462 fFWExtension = 0;
464}
465
466////////////////////////////////////////////////////////////////////////////////
467/// default constructor
468
469TGeoVolume::TGeoVolume(const char *name, const TGeoShape *shape, const TGeoMedium *med)
470 :TNamed(name, "")
471{
472 fName = fName.Strip();
473 fNodes = 0;
474 fShape = (TGeoShape*)shape;
475 if (fShape) {
477 Warning("Ctor", "volume %s has invalid shape", name);
478 }
479 if (!fShape->IsValid()) {
480 Fatal("ctor", "Shape of volume %s invalid. Aborting!", fName.Data());
481 }
482 }
483 fMedium = (TGeoMedium*)med;
485 fFinder = 0;
486 fVoxels = 0;
488 fField = 0;
489 fOption = "";
490 fNumber = 0;
491 fNtotal = 0;
492 fRefCount = 0;
493 fUserExtension = 0;
494 fFWExtension = 0;
497}
498
499////////////////////////////////////////////////////////////////////////////////
500/// Destructor
501
503{
504 if (fNodes) {
506 fNodes->Delete();
507 }
508 delete fNodes;
509 }
511 if (fVoxels) delete fVoxels;
514}
515
516////////////////////////////////////////////////////////////////////////////////
517/// How to browse a volume
518
520{
521 if (!b) return;
522
523// if (!GetNdaughters()) b->Add(this, GetName(), IsVisible());
524 TGeoVolume *daughter;
525 TString title;
526 for (Int_t i=0; i<GetNdaughters(); i++) {
527 daughter = GetNode(i)->GetVolume();
528 if(daughter->GetTitle()[0]) {
529 if (daughter->IsAssembly()) title.TString::Format("Assembly with %d daughter(s)",
530 daughter->GetNdaughters());
531 else if (daughter->GetFinder()) {
532 TString s1 = daughter->GetFinder()->ClassName();
533 s1.ReplaceAll("TGeoPattern","");
534 title.TString::Format("Volume having %s shape divided in %d %s slices",
535 daughter->GetShape()->ClassName(),daughter->GetNdaughters(), s1.Data());
536
537 } else title.TString::Format("Volume with %s shape having %d daughter(s)",
538 daughter->GetShape()->ClassName(),daughter->GetNdaughters());
539 daughter->SetTitle(title.Data());
540 }
541 b->Add(daughter, daughter->GetName(), daughter->IsVisible());
542// if (IsVisDaughters())
543// b->AddCheckBox(daughter, daughter->IsVisible());
544// else
545// b->AddCheckBox(daughter, kFALSE);
546 }
547}
548
549////////////////////////////////////////////////////////////////////////////////
550/// Computes the capacity of this [cm^3] as the capacity of its shape.
551/// In case of assemblies, the capacity is computed as the sum of daughter's capacities.
552
554{
555 if (!IsAssembly()) return fShape->Capacity();
556 Double_t capacity = 0.0;
557 Int_t nd = GetNdaughters();
558 Int_t i;
559 for (i=0; i<nd; i++) capacity += GetNode(i)->GetVolume()->Capacity();
560 return capacity;
561}
562
563////////////////////////////////////////////////////////////////////////////////
564/// Shoot nrays with random directions from starting point (startx, starty, startz)
565/// in the reference frame of this volume. Track each ray until exiting geometry, then
566/// shoot backwards from exiting point and compare boundary crossing points.
567
568void TGeoVolume::CheckGeometry(Int_t nrays, Double_t startx, Double_t starty, Double_t startz) const
569{
570 TGeoVolume *old_vol = fGeoManager->GetTopVolume();
571 if (old_vol!=this) fGeoManager->SetTopVolume((TGeoVolume*)this);
572 else old_vol=0;
575 painter->CheckGeometry(nrays, startx, starty, startz);
576}
577
578////////////////////////////////////////////////////////////////////////////////
579/// Overlap checking tool. Check for illegal overlaps within a limit OVLP.
580/// Use option="s[number]" to force overlap checking by sampling volume with
581/// [number] points.
582///
583/// Ex:
584/// ~~~ {.cpp}
585/// myVol->CheckOverlaps(0.01, "s10000000"); // shoot 10000000 points
586/// myVol->CheckOverlaps(0.01, "s"); // shoot the default value of 1e6 points
587/// ~~~
588
590{
591 if (!GetNdaughters() || fFinder) return;
592 Bool_t sampling = kFALSE;
593 TString opt(option);
594 opt.ToLower();
595 if (opt.Contains("s")) sampling = kTRUE;
597 if (!sampling) fGeoManager->SetNsegments(80);
600// Info("CheckOverlaps", "=== Checking overlaps for volume %s ===\n", GetName());
601 }
602 painter->CheckOverlaps(this, ovlp, option);
603// if (sampling) return;
607 Int_t novlps = overlaps->GetEntriesFast();
608 TNamed *obj;
610 for (Int_t i=0; i<novlps; i++) {
611 obj = (TNamed*)overlaps->At(i);
612 if (novlps<1000) name = TString::Format("ov%03d", i);
613 else name = TString::Format("ov%06d", i);
614 obj->SetName(name);
615 }
616 if (novlps) Info("CheckOverlaps", "Number of illegal overlaps/extrusions for volume %s: %d\n", GetName(), novlps);
617 }
618}
619
620////////////////////////////////////////////////////////////////////////////////
621/// Tests for checking the shape navigation algorithms. See TGeoShape::CheckShape()
622
624{
625 fShape->CheckShape(testNo,nsamples,option);
626}
627
628////////////////////////////////////////////////////////////////////////////////
629/// Clean data of the volume.
630
632{
633 ClearNodes();
634 ClearShape();
635}
636
637////////////////////////////////////////////////////////////////////////////////
638/// Clear the shape of this volume from the list held by the current manager.
639
641{
643}
644
645////////////////////////////////////////////////////////////////////////////////
646/// check for negative parameters in shapes.
647
649{
650 if (fShape->IsRunTimeShape()) {
651 Error("CheckShapes", "volume %s has run-time shape", GetName());
652 InspectShape();
653 return;
654 }
655 if (!fNodes) return;
657 TGeoNode *node = 0;
658 TGeoNode *new_node;
659 const TGeoShape *shape = 0;
660 TGeoVolume *old_vol;
661 for (Int_t i=0; i<nd; i++) {
662 node=(TGeoNode*)fNodes->At(i);
663 // check if node has name
664 if (!node->GetName()[0]) printf("Daughter %i of volume %s - NO NAME!!!\n",
665 i, GetName());
666 old_vol = node->GetVolume();
667 shape = old_vol->GetShape();
668 if (shape->IsRunTimeShape()) {
669// printf(" Node %s/%s has shape with negative parameters. \n",
670// GetName(), node->GetName());
671// old_vol->InspectShape();
672 // make a copy of the node
673 new_node = node->MakeCopyNode();
674 if (!new_node) {
675 Fatal("CheckShapes", "Cannot make copy node for %s", node->GetName());
676 return;
677 }
678 TGeoShape *new_shape = shape->GetMakeRuntimeShape(fShape, node->GetMatrix());
679 if (!new_shape) {
680 Error("CheckShapes","cannot resolve runtime shape for volume %s/%s\n",
681 GetName(),old_vol->GetName());
682 continue;
683 }
684 TGeoVolume *new_volume = old_vol->MakeCopyVolume(new_shape);
685// printf(" new volume %s shape params :\n", new_volume->GetName());
686// new_volume->InspectShape();
687 new_node->SetVolume(new_volume);
688 // decouple the old node and put the new one instead
689 fNodes->AddAt(new_node, i);
690// new_volume->CheckShapes();
691 }
692 }
693}
694
695////////////////////////////////////////////////////////////////////////////////
696/// Count total number of subnodes starting from this volume, nlevels down
697/// - option = 0 (default) - count only once per volume
698/// - option = 1 - count every time
699/// - option = 2 - count volumes on visible branches
700/// - option = 3 - return maximum level counted already with option = 0
701
703{
704 static Int_t maxlevel = 0;
705 static Int_t nlev = 0;
706
707 if (option<0 || option>3) option = 0;
708 Int_t visopt = 0;
709 Int_t nd = GetNdaughters();
710 Bool_t last = (!nlevels || !nd)?kTRUE:kFALSE;
711 switch (option) {
712 case 0:
713 if (fNtotal) return fNtotal;
714 case 1:
715 fNtotal = 1;
716 break;
717 case 2:
718 visopt = fGeoManager->GetVisOption();
719 if (!IsVisDaughters()) last = kTRUE;
720 switch (visopt) {
722 fNtotal = (IsVisible())?1:0;
723 break;
725 fNtotal = (IsVisible() && last)?1:0;
726 }
727 if (!IsVisibleDaughters()) return fNtotal;
728 break;
729 case 3:
730 return maxlevel;
731 }
732 if (last) return fNtotal;
733 if (gGeoManager->GetTopVolume() == this) {
734 maxlevel=0;
735 nlev = 0;
736 }
737 if (nlev>maxlevel) maxlevel = nlev;
738 TGeoNode *node;
739 TGeoVolume *vol;
740 nlev++;
741 for (Int_t i=0; i<nd; i++) {
742 node = GetNode(i);
743 vol = node->GetVolume();
744 fNtotal += vol->CountNodes(nlevels-1, option);
745 }
746 nlev--;
747 return fNtotal;
748}
749
750////////////////////////////////////////////////////////////////////////////////
751/// Return TRUE if volume and all daughters are invisible.
752
754{
755 if (IsVisible()) return kFALSE;
756 Int_t nd = GetNdaughters();
757 for (Int_t i=0; i<nd; i++) if (GetNode(i)->GetVolume()->IsVisible()) return kFALSE;
758 return kTRUE;
759}
760
761////////////////////////////////////////////////////////////////////////////////
762/// Make volume and each of it daughters (in)visible.
763
765{
766 SetAttVisibility(!flag);
767 Int_t nd = GetNdaughters();
768 TObjArray *list = new TObjArray(nd+1);
769 list->Add(this);
770 TGeoVolume *vol;
771 for (Int_t i=0; i<nd; i++) {
772 vol = GetNode(i)->GetVolume();
773 vol->SetAttVisibility(!flag);
774 list->Add(vol);
775 }
776 TIter next(gROOT->GetListOfBrowsers());
777 TBrowser *browser = 0;
778 while ((browser=(TBrowser*)next())) {
779 for (Int_t i=0; i<nd+1; i++) {
780 vol = (TGeoVolume*)list->At(i);
781 browser->CheckObjectItem(vol, !flag);
782 }
783 browser->Refresh();
784 }
785 delete list;
787}
788
789////////////////////////////////////////////////////////////////////////////////
790/// Return TRUE if volume contains nodes
791
793{
794 return kTRUE;
795}
796
797////////////////////////////////////////////////////////////////////////////////
798/// check if the visibility and attributes are the default ones
799
801{
802 if (!IsVisible()) return kFALSE;
803 if (GetLineColor() != gStyle->GetLineColor()) return kFALSE;
804 if (GetLineStyle() != gStyle->GetLineStyle()) return kFALSE;
805 if (GetLineWidth() != gStyle->GetLineWidth()) return kFALSE;
806 return kTRUE;
807}
808
809////////////////////////////////////////////////////////////////////////////////
810/// True if this is the top volume of the geometry
811
813{
814 if (fGeoManager->GetTopVolume() == this) return kTRUE;
815 return kFALSE;
816}
817
818////////////////////////////////////////////////////////////////////////////////
819/// Check if the painter is currently ray-tracing the content of this volume.
820
822{
823 return TGeoAtt::IsVisRaytrace();
824}
825
826////////////////////////////////////////////////////////////////////////////////
827/// Inspect the material for this volume.
828
830{
831 GetMaterial()->Print();
832}
833
834////////////////////////////////////////////////////////////////////////////////
835/// Import a volume from a file.
836
837TGeoVolume *TGeoVolume::Import(const char *filename, const char *name, Option_t * /*option*/)
838{
839 if (!gGeoManager) gGeoManager = new TGeoManager("geometry","");
840 if (!filename) return 0;
841 TGeoVolume *volume = 0;
842 if (strstr(filename,".gdml")) {
843 // import from a gdml file
844 } else {
845 // import from a root file
848 if (!f || f->IsZombie()) {
849 printf("Error: TGeoVolume::Import : Cannot open file %s\n", filename);
850 return 0;
851 }
852 if (name && name[0]) {
853 volume = (TGeoVolume*)f->Get(name);
854 } else {
855 TIter next(f->GetListOfKeys());
856 TKey *key;
857 while ((key = (TKey*)next())) {
858 if (strcmp(key->GetClassName(),"TGeoVolume") != 0) continue;
859 volume = (TGeoVolume*)key->ReadObj();
860 break;
861 }
862 }
863 delete f;
864 }
865 if (!volume) return NULL;
866 volume->RegisterYourself();
867 return volume;
868}
869
870////////////////////////////////////////////////////////////////////////////////
871/// Export this volume to a file.
872///
873/// - Case 1: root file or root/xml file
874/// if filename end with ".root". The key will be named name
875/// if filename end with ".xml" a root/xml file is produced.
876///
877/// - Case 2: C++ script
878/// if filename end with ".C"
879///
880/// - Case 3: gdml file
881/// if filename end with ".gdml"
882///
883/// NOTE that to use this option, the PYTHONPATH must be defined like
884/// export PYTHONPATH=$ROOTSYS/lib:$ROOTSYS/gdml
885///
886
888{
889 TString sfile(filename);
890 if (sfile.Contains(".C")) {
891 //Save volume as a C++ script
892 Info("Export","Exporting volume %s as C++ code", GetName());
893 SaveAs(filename, "");
894 return 1;
895 }
896 if (sfile.Contains(".gdml")) {
897 //Save geometry as a gdml file
898 Info("Export","Exporting %s as gdml code - not implemented yet", GetName());
899 return 0;
900 }
901 if (sfile.Contains(".root") || sfile.Contains(".xml")) {
902 //Save volume in a root file
903 Info("Export","Exporting %s as root file.", GetName());
904 TString opt(option);
905 if (!opt.Length()) opt = "recreate";
906 TFile *f = TFile::Open(filename,opt.Data());
907 if (!f || f->IsZombie()) {
908 Error("Export","Cannot open file");
909 return 0;
910 }
911 TString keyname(name);
912 if (keyname.IsNull()) keyname = GetName();
913 Int_t nbytes = Write(keyname);
914 delete f;
915 return nbytes;
916 }
917 return 0;
918}
919
920////////////////////////////////////////////////////////////////////////////////
921/// Actualize matrix of node indexed `<inode>`
922
923void TGeoVolume::cd(Int_t inode) const
924{
925 if (fFinder) fFinder->cd(inode-fFinder->GetDivIndex());
926}
927
928////////////////////////////////////////////////////////////////////////////////
929/// Add a TGeoNode to the list of nodes. This is the usual method for adding
930/// daughters inside the container volume.
931
933{
934 TGeoMatrix *matrix = mat;
935 if (matrix==0) matrix = gGeoIdentity;
936 else matrix->RegisterYourself();
937 if (!vol) {
938 Error("AddNode", "Volume is NULL");
939 return 0;
940 }
941 if (!vol->IsValid()) {
942 Error("AddNode", "Won't add node with invalid shape");
943 printf("### invalid volume was : %s\n", vol->GetName());
944 return 0;
945 }
946 if (!fNodes) fNodes = new TObjArray();
947
948 if (fFinder) {
949 // volume already divided.
950 Error("AddNode", "Cannot add node %s_%i into divided volume %s", vol->GetName(), copy_no, GetName());
951 return 0;
952 }
953
954 TGeoNodeMatrix *node = 0;
955 node = new TGeoNodeMatrix(vol, matrix);
956 node->SetMotherVolume(this);
957 fNodes->Add(node);
958 TString name = TString::Format("%s_%d", vol->GetName(), copy_no);
959// if (fNodes->FindObject(name))
960// Warning("AddNode", "Volume %s : added node %s with same name", GetName(), name.Data());
961 node->SetName(name);
962 node->SetNumber(copy_no);
963 fRefCount++;
964 vol->Grab();
965 return node;
966}
967
968////////////////////////////////////////////////////////////////////////////////
969/// Add a division node to the list of nodes. The method is called by
970/// TGeoVolume::Divide() for creating the division nodes.
971
973{
974 if (!vol) {
975 Error("AddNodeOffset", "invalid volume");
976 return;
977 }
978 if (!vol->IsValid()) {
979 Error("AddNode", "Won't add node with invalid shape");
980 printf("### invalid volume was : %s\n", vol->GetName());
981 return;
982 }
983 if (!fNodes) fNodes = new TObjArray();
984 TGeoNode *node = new TGeoNodeOffset(vol, copy_no, offset);
985 node->SetMotherVolume(this);
986 fNodes->Add(node);
987 TString name = TString::Format("%s_%d", vol->GetName(), copy_no+1);
988 node->SetName(name);
989 node->SetNumber(copy_no+1);
990 vol->Grab();
991}
992
993////////////////////////////////////////////////////////////////////////////////
994/// Add a TGeoNode to the list of nodes. This is the usual method for adding
995/// daughters inside the container volume.
996
998{
999 if (!vol) {
1000 Error("AddNodeOverlap", "Volume is NULL");
1001 return;
1002 }
1003 if (!vol->IsValid()) {
1004 Error("AddNodeOverlap", "Won't add node with invalid shape");
1005 printf("### invalid volume was : %s\n", vol->GetName());
1006 return;
1007 }
1008 if (vol->IsAssembly()) {
1009 Warning("AddNodeOverlap", "Declaring assembly %s as possibly overlapping inside %s not allowed. Using AddNode instead !",vol->GetName(),GetName());
1010 AddNode(vol, copy_no, mat, option);
1011 return;
1012 }
1013 TGeoMatrix *matrix = mat;
1014 if (matrix==0) matrix = gGeoIdentity;
1015 else matrix->RegisterYourself();
1016 if (!fNodes) fNodes = new TObjArray();
1017
1018 if (fFinder) {
1019 // volume already divided.
1020 Error("AddNodeOverlap", "Cannot add node %s_%i into divided volume %s", vol->GetName(), copy_no, GetName());
1021 return;
1022 }
1023
1024 TGeoNodeMatrix *node = new TGeoNodeMatrix(vol, matrix);
1025 node->SetMotherVolume(this);
1026 fNodes->Add(node);
1027 TString name = TString::Format("%s_%d", vol->GetName(), copy_no);
1028 if (fNodes->FindObject(name))
1029 Warning("AddNode", "Volume %s : added node %s with same name", GetName(), name.Data());
1030 node->SetName(name);
1031 node->SetNumber(copy_no);
1032 node->SetOverlapping();
1033 if (vol->GetMedium() == fMedium)
1034 node->SetVirtual();
1035 vol->Grab();
1036}
1037
1038////////////////////////////////////////////////////////////////////////////////
1039/// Division a la G3. The volume will be divided along IAXIS (see shape classes), in NDIV
1040/// slices, from START with given STEP. The division volumes will have medium number NUMED.
1041/// If NUMED=0 they will get the medium number of the divided volume (this). If NDIV<=0,
1042/// all range of IAXIS will be divided and the resulting number of divisions will be centered on
1043/// IAXIS. If STEP<=0, the real STEP will be computed as the full range of IAXIS divided by NDIV.
1044/// Options (case insensitive):
1045/// - N - divide all range in NDIV cells (same effect as STEP<=0) (GSDVN in G3)
1046/// - NX - divide range starting with START in NDIV cells (GSDVN2 in G3)
1047/// - S - divide all range with given STEP. NDIV is computed and divisions will be centered
1048/// in full range (same effect as NDIV<=0) (GSDVS, GSDVT in G3)
1049/// - SX - same as DVS, but from START position. (GSDVS2, GSDVT2 in G3)
1050
1051TGeoVolume *TGeoVolume::Divide(const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step, Int_t numed, Option_t *option)
1052{
1053 if (fFinder) {
1054 // volume already divided.
1055 Fatal("Divide","volume %s already divided", GetName());
1056 return 0;
1057 }
1058 TString opt(option);
1059 opt.ToLower();
1060 TString stype = fShape->ClassName();
1061 if (!fNodes) fNodes = new TObjArray();
1062 Double_t xlo, xhi, range;
1063 range = fShape->GetAxisRange(iaxis, xlo, xhi);
1064 // for phi divisions correct the range
1065 if (!strcmp(fShape->GetAxisName(iaxis), "PHI")) {
1066 if ((start-xlo)<-1E-3) start+=360.;
1067 if (TGeoShape::IsSameWithinTolerance(range,360)) {
1068 xlo = start;
1069 xhi = start+range;
1070 }
1071 }
1072 if (range <=0) {
1073 InspectShape();
1074 Fatal("Divide", "cannot divide volume %s (%s) on %s axis", GetName(), stype.Data(), fShape->GetAxisName(iaxis));
1075 return 0;
1076 }
1077 if (ndiv<=0 || opt.Contains("s")) {
1078 if (step<=0) {
1079 Fatal("Divide", "invalid division type for volume %s : ndiv=%i, step=%g", GetName(), ndiv, step);
1080 return 0;
1081 }
1082 if (opt.Contains("x")) {
1083 if ((xlo-start)>1E-3 || (xhi-start)<-1E-3) {
1084 Fatal("Divide", "invalid START=%g for division on axis %s of volume %s. Range is (%g, %g)",
1085 start, fShape->GetAxisName(iaxis), GetName(), xlo, xhi);
1086 return 0;
1087 }
1088 xlo = start;
1089 range = xhi-xlo;
1090 }
1091 ndiv = Int_t((range+0.1*step)/step);
1092 Double_t ddx = range - ndiv*step;
1093 // always center the division in this case
1094 if (ddx>1E-3) Warning("Divide", "division of volume %s on %s axis (ndiv=%d) will be centered in the full range",
1095 GetName(), fShape->GetAxisName(iaxis), ndiv);
1096 start = xlo + 0.5*ddx;
1097 }
1098 if (step<=0 || opt.Contains("n")) {
1099 if (opt.Contains("x")) {
1100 if ((xlo-start)>1E-3 || (xhi-start)<-1E-3) {
1101 Fatal("Divide", "invalid START=%g for division on axis %s of volume %s. Range is (%g, %g)",
1102 start, fShape->GetAxisName(iaxis), GetName(), xlo, xhi);
1103 return 0;
1104 }
1105 xlo = start;
1106 range = xhi-xlo;
1107 }
1108 step = range/ndiv;
1109 start = xlo;
1110 }
1111
1112 Double_t end = start+ndiv*step;
1113 if (((start-xlo)<-1E-3) || ((end-xhi)>1E-3)) {
1114 Fatal("Divide", "division of volume %s on axis %s exceed range (%g, %g)",
1115 GetName(), fShape->GetAxisName(iaxis), xlo, xhi);
1116 return 0;
1117 }
1118 TGeoVolume *voldiv = fShape->Divide(this, divname, iaxis, ndiv, start, step);
1119 if (numed) {
1120 TGeoMedium *medium = fGeoManager->GetMedium(numed);
1121 if (!medium) {
1122 Fatal("Divide", "invalid medium number %d for division volume %s", numed, divname);
1123 return voldiv;
1124 }
1125 voldiv->SetMedium(medium);
1126 if (medium->GetMaterial()) medium->GetMaterial()->SetUsed();
1127 }
1128 return voldiv;
1129}
1130
1131////////////////////////////////////////////////////////////////////////////////
1132/// compute the closest distance of approach from point px,py to this volume
1133
1135{
1138 Int_t dist = 9999;
1139 if (!painter) return dist;
1140 dist = painter->DistanceToPrimitiveVol(this, px, py);
1141 return dist;
1142}
1143
1144////////////////////////////////////////////////////////////////////////////////
1145/// draw top volume according to option
1146
1148{
1153 if (!IsVisContainers()) SetVisLeaves();
1154 if (option && option[0] > 0) {
1155 painter->DrawVolume(this, option);
1156 } else {
1157 painter->DrawVolume(this, gEnv->GetValue("Viewer3D.DefaultDrawOption",""));
1158 }
1159}
1160
1161////////////////////////////////////////////////////////////////////////////////
1162/// draw only this volume
1163
1165{
1166 if (IsAssembly()) {
1167 Info("DrawOnly", "Volume assemblies do not support this option.");
1168 return;
1169 }
1171 SetVisOnly();
1174 if (option && option[0] > 0) {
1175 painter->DrawVolume(this, option);
1176 } else {
1177 painter->DrawVolume(this, gEnv->GetValue("Viewer3D.DefaultDrawOption",""));
1178 }
1179}
1180
1181////////////////////////////////////////////////////////////////////////////////
1182/// Perform an extensive sampling to find which type of voxelization is
1183/// most efficient.
1184
1186{
1187 printf("Optimizing volume %s ...\n", GetName());
1189 return painter->TestVoxels(this);
1190}
1191
1192////////////////////////////////////////////////////////////////////////////////
1193/// Print volume info
1194
1196{
1197 printf("== Volume: %s type %s positioned %d times\n", GetName(), ClassName(), fRefCount);
1198 InspectShape();
1200}
1201
1202////////////////////////////////////////////////////////////////////////////////
1203/// paint volume
1204
1206{
1208 painter->SetTopVolume(this);
1209// painter->Paint(option);
1210 if (option && option[0] > 0) {
1211 painter->Paint(option);
1212 } else {
1213 painter->Paint(gEnv->GetValue("Viewer3D.DefaultDrawOption",""));
1214 }
1215}
1216
1217////////////////////////////////////////////////////////////////////////////////
1218/// Print the voxels for this volume.
1219
1221{
1222 if (fVoxels) fVoxels->Print();
1223}
1224
1225////////////////////////////////////////////////////////////////////////////////
1226/// Recreate the content of the other volume without pointer copying. Voxels are
1227/// ignored and supposed to be created in a later step via Voxelize.
1228
1230{
1231 Int_t nd = other->GetNdaughters();
1232 if (!nd) return;
1233 TGeoPatternFinder *finder = other->GetFinder();
1234 if (finder) {
1235 Int_t iaxis = finder->GetDivAxis();
1236 Int_t ndiv = finder->GetNdiv();
1237 Double_t start = finder->GetStart();
1238 Double_t step = finder->GetStep();
1239 Int_t numed = other->GetNode(0)->GetVolume()->GetMedium()->GetId();
1240 TGeoVolume *voldiv = Divide(other->GetNode(0)->GetVolume()->GetName(), iaxis, ndiv, start, step, numed);
1241 voldiv->ReplayCreation(other->GetNode(0)->GetVolume());
1242 return;
1243 }
1244 for (Int_t i=0; i<nd; i++) {
1245 TGeoNode *node = other->GetNode(i);
1246 if (node->IsOverlapping()) AddNodeOverlap(node->GetVolume(), node->GetNumber(), node->GetMatrix());
1247 else AddNode(node->GetVolume(), node->GetNumber(), node->GetMatrix());
1248 }
1249}
1250
1251////////////////////////////////////////////////////////////////////////////////
1252/// print nodes
1253
1255{
1256 Int_t nd = GetNdaughters();
1257 for (Int_t i=0; i<nd; i++) {
1258 printf("%s\n", GetNode(i)->GetName());
1259 cd(i);
1260 GetNode(i)->GetMatrix()->Print();
1261 }
1262}
1263////////////////////////////////////////////////////////////////////////////////
1264/// Generate a lego plot fot the top volume, according to option.
1265
1267 Int_t nphi, Double_t phimin, Double_t phimax,
1268 Double_t rmin, Double_t rmax, Option_t *option)
1269{
1271 TGeoVolume *old_vol = fGeoManager->GetTopVolume();
1272 if (old_vol!=this) fGeoManager->SetTopVolume(this);
1273 else old_vol=0;
1274 TH2F *hist = p->LegoPlot(ntheta, themin, themax, nphi, phimin, phimax, rmin, rmax, option);
1275 hist->Draw("lego1sph");
1276 return hist;
1277}
1278
1279////////////////////////////////////////////////////////////////////////////////
1280/// Register the volume and all materials/media/matrices/shapes to the manager.
1281
1283{
1284 if (fGeoManager->GetListOfVolumes()->FindObject(this)) return;
1285 // Register volume
1286 fGeoManager->AddVolume(this);
1287 // Register shape
1289 if (fShape->IsComposite()) {
1291 comp->RegisterYourself();
1292 } else {
1294 }
1295 }
1296 // Register medium/material
1301 }
1302 // Register matrices for nodes.
1303 TGeoMatrix *matrix;
1304 TGeoNode *node;
1305 Int_t nd = GetNdaughters();
1306 Int_t i;
1307 for (i=0; i<nd; i++) {
1308 node = GetNode(i);
1309 matrix = node->GetMatrix();
1310 if (!matrix->IsRegistered()) matrix->RegisterYourself();
1311 else if (!fGeoManager->GetListOfMatrices()->FindObject(matrix)) {
1312 fGeoManager->GetListOfMatrices()->Add(matrix);
1313 }
1314 }
1315 // Call RegisterYourself recursively
1316 for (i=0; i<nd; i++) GetNode(i)->GetVolume()->RegisterYourself(option);
1317}
1318
1319////////////////////////////////////////////////////////////////////////////////
1320/// Draw random points in the bounding box of this volume.
1321
1323{
1325 TGeoVolume *old_vol = fGeoManager->GetTopVolume();
1326 if (old_vol!=this) fGeoManager->SetTopVolume(this);
1327 else old_vol=0;
1328 fGeoManager->RandomPoints(this, npoints, option);
1329 if (old_vol) fGeoManager->SetTopVolume(old_vol);
1330}
1331
1332////////////////////////////////////////////////////////////////////////////////
1333/// Random raytracing method.
1334
1335void TGeoVolume::RandomRays(Int_t nrays, Double_t startx, Double_t starty, Double_t startz, const char *target_vol, Bool_t check_norm)
1336{
1338 TGeoVolume *old_vol = fGeoManager->GetTopVolume();
1339 if (old_vol!=this) fGeoManager->SetTopVolume(this);
1340 else old_vol=0;
1341 fGeoManager->RandomRays(nrays, startx, starty, startz, target_vol, check_norm);
1342 if (old_vol) fGeoManager->SetTopVolume(old_vol);
1343}
1344
1345////////////////////////////////////////////////////////////////////////////////
1346/// Draw this volume with current settings and perform raytracing in the pad.
1347
1349{
1353 Bool_t drawn = (painter->GetDrawnVolume()==this)?kTRUE:kFALSE;
1354 if (!drawn) {
1355 painter->DrawVolume(this, "");
1357 painter->ModifiedPad();
1358 return;
1359 }
1361 painter->ModifiedPad();
1362}
1363
1364////////////////////////////////////////////////////////////////////////////////
1365/// Save geometry having this as top volume as a C++ macro.
1366
1368{
1369 if (!filename) return;
1370 std::ofstream out;
1371 out.open(filename, std::ios::out);
1372 if (out.bad()) {
1373 Error("SavePrimitive", "Bad file name: %s", filename);
1374 return;
1375 }
1377
1378 TString fname(filename);
1379 Int_t ind = fname.Index(".");
1380 if (ind>0) fname.Remove(ind);
1381 out << "void "<<fname<<"() {" << std::endl;
1382 out << " gSystem->Load(\"libGeom\");" << std::endl;
1384 out << std::setprecision(prec);
1385 ((TGeoVolume*)this)->SavePrimitive(out,option);
1386 out << "}" << std::endl;
1387}
1388
1389////////////////////////////////////////////////////////////////////////////////
1390/// Connect user-defined extension to the volume. The volume "grabs" a copy, so
1391/// the original object can be released by the producer. Release the previously
1392/// connected extension if any.
1393///
1394/// NOTE: This interface is intended for user extensions and is guaranteed not
1395/// to be used by TGeo
1396
1398{
1400 fUserExtension = 0;
1401 if (ext) fUserExtension = ext->Grab();
1402}
1403
1404////////////////////////////////////////////////////////////////////////////////
1405/// Connect framework defined extension to the volume. The volume "grabs" a copy,
1406/// so the original object can be released by the producer. Release the previously
1407/// connected extension if any.
1408///
1409/// NOTE: This interface is intended for the use by TGeo and the users should
1410/// NOT connect extensions using this method
1411
1413{
1415 fFWExtension = 0;
1416 if (ext) fFWExtension = ext->Grab();
1417}
1418
1419////////////////////////////////////////////////////////////////////////////////
1420/// Get a copy of the user extension pointer. The user must call Release() on
1421/// the copy pointer once this pointer is not needed anymore (equivalent to
1422/// delete() after calling new())
1423
1425{
1426 if (fUserExtension) return fUserExtension->Grab();
1427 return 0;
1428}
1429
1430////////////////////////////////////////////////////////////////////////////////
1431/// Get a copy of the framework extension pointer. The user must call Release() on
1432/// the copy pointer once this pointer is not needed anymore (equivalent to
1433/// delete() after calling new())
1434
1436{
1437 if (fFWExtension) return fFWExtension->Grab();
1438 return 0;
1439}
1440
1441////////////////////////////////////////////////////////////////////////////////
1442/// Save a primitive as a C++ statement(s) on output stream "out".
1443
1444void TGeoVolume::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
1445{
1446 Int_t i,icopy;
1447 Int_t nd = GetNdaughters();
1448 TGeoVolume *dvol;
1449 TGeoNode *dnode;
1450 TGeoMatrix *matrix;
1451
1452 // check if we need to save shape/volume
1453 Bool_t mustDraw = kFALSE;
1454 if (fGeoManager->GetGeomPainter()->GetTopVolume()==this) mustDraw = kTRUE;
1455 if (!option[0]) {
1457 out << " new TGeoManager(\"" << fGeoManager->GetName() << "\", \"" << fGeoManager->GetTitle() << "\");" << std::endl << std::endl;
1458// if (mustDraw) out << " Bool_t mustDraw = kTRUE;" << std::endl;
1459// else out << " Bool_t mustDraw = kFALSE;" << std::endl;
1460 out << " Double_t dx, dy, dz;" << std::endl;
1461 out << " Double_t dx1, dx2, dy1, dy2;" << std::endl;
1462 out << " Double_t vert[20], par[20];" << std::endl;
1463 out << " Double_t theta, phi, h1, bl1, tl1, alpha1, h2, bl2, tl2, alpha2;" << std::endl;
1464 out << " Double_t twist;" << std::endl;
1465 out << " Double_t origin[3];" << std::endl;
1466 out << " Double_t rmin, rmax, rmin1, rmax1, rmin2, rmax2;" << std::endl;
1467 out << " Double_t r, rlo, rhi;" << std::endl;
1468 out << " Double_t a, b;" << std::endl;
1469 out << " Double_t point[3], norm[3];" << std::endl;
1470 out << " Double_t rin, stin, rout, stout;" << std::endl;
1471 out << " Double_t thx, phx, thy, phy, thz, phz;" << std::endl;
1472 out << " Double_t alpha, theta1, theta2, phi1, phi2, dphi;" << std::endl;
1473 out << " Double_t tr[3], rot[9];" << std::endl;
1474 out << " Double_t z, density, radl, absl, w;" << std::endl;
1475 out << " Double_t lx, ly, lz, tx, ty, tz;" << std::endl;
1476 out << " Double_t xvert[50], yvert[50];" << std::endl;
1477 out << " Double_t zsect, x0, y0, scale0;" << std::endl;
1478 out << " Int_t nel, numed, nz, nedges, nvert;" << std::endl;
1479 out << " TGeoBoolNode *pBoolNode = nullptr;" << std::endl << std::endl;
1480 // first save materials/media
1481 out << " // MATERIALS, MIXTURES AND TRACKING MEDIA" << std::endl;
1482 SavePrimitive(out, "m");
1483 // then, save matrices
1484 out << std::endl << " // TRANSFORMATION MATRICES" << std::endl;
1485 SavePrimitive(out, "x");
1486 // save this volume and shape
1487 SavePrimitive(out, "s");
1488 out << std::endl << " // SET TOP VOLUME OF GEOMETRY" << std::endl;
1489 out << " gGeoManager->SetTopVolume(" << GetPointerName() << ");" << std::endl;
1490 // save daughters
1491 out << std::endl << " // SHAPES, VOLUMES AND GEOMETRICAL HIERARCHY" << std::endl;
1492 SavePrimitive(out, "d");
1493 out << std::endl << " // CLOSE GEOMETRY" << std::endl;
1494 out << " gGeoManager->CloseGeometry();" << std::endl;
1495 if (mustDraw) {
1496 if (!IsRaytracing()) out << " gGeoManager->GetTopVolume()->Draw();" << std::endl;
1497 else out << " gGeoManager->GetTopVolume()->Raytrace();" << std::endl;
1498 }
1499 return;
1500 }
1501 // check if we need to save shape/volume
1502 if (!strcmp(option, "s")) {
1503 // create the shape for this volume
1505 if (!IsAssembly()) {
1507 out << " // Volume: " << GetName() << std::endl;
1508 out << " TGeoVolume *" << GetPointerName() << " = new TGeoVolume(\"" << GetName() << "\"," << fShape->GetPointerName();
1509 if (fMedium) out << ", " << fMedium->GetPointerName();
1510 out << ");" << std::endl;
1511 } else {
1512 out << " // Assembly: " << GetName() << std::endl;
1513 out << " " << GetPointerName() << " = new TGeoVolumeAssembly(\"" << GetName() << "\"" << ");" << std::endl;
1514 }
1515 SaveLineAttributes(out, GetPointerName(), 1, 1, 1);
1516 if (!IsVisible() && !IsAssembly()) out << " " << GetPointerName() << "->SetVisibility(kFALSE);" << std::endl;
1517 if (!IsVisibleDaughters()) out << " " << GetPointerName() << "->VisibleDaughters(kFALSE);" << std::endl;
1518 if (IsVisContainers()) out << " " << GetPointerName() << "->SetVisContainers(kTRUE);" << std::endl;
1519 if (IsVisLeaves()) out << " " << GetPointerName() << "->SetVisLeaves(kTRUE);" << std::endl;
1521 }
1522 // check if we need to save the media
1523 if (!strcmp(option, "m")) {
1525 for (i=0; i<nd; i++) {
1526 dvol = GetNode(i)->GetVolume();
1527 dvol->SavePrimitive(out,option);
1528 }
1529 return;
1530 }
1531 // check if we need to save the matrices
1532 if (!strcmp(option, "x")) {
1533 if (fFinder) {
1534 dvol = GetNode(0)->GetVolume();
1535 dvol->SavePrimitive(out,option);
1536 return;
1537 }
1538 for (i=0; i<nd; i++) {
1539 dnode = GetNode(i);
1540 matrix = dnode->GetMatrix();
1541 if (!matrix->IsIdentity()) matrix->SavePrimitive(out,option);
1542 dnode->GetVolume()->SavePrimitive(out,option);
1543 }
1544 return;
1545 }
1546 // check if we need to save volume daughters
1547 if (!strcmp(option, "d")) {
1548 if (!nd) return;
1551 if (fFinder) {
1552 // volume divided: generate volume->Divide()
1553 dnode = GetNode(0);
1554 dvol = dnode->GetVolume();
1555 out << " TGeoVolume *" << dvol->GetPointerName() << " = ";
1556 out << GetPointerName() << "->Divide(\"" << dvol->GetName() << "\", ";
1558 if (fMedium != dvol->GetMedium())
1559 out << ", " << dvol->GetMedium()->GetId();
1560 out << ");" << std::endl;
1561 dvol->SavePrimitive(out,"d");
1562 return;
1563 }
1564 for (i=0; i<nd; i++) {
1565 dnode = GetNode(i);
1566 dvol = dnode->GetVolume();
1567 dvol->SavePrimitive(out,"s");
1568 matrix = dnode->GetMatrix();
1569 icopy = dnode->GetNumber();
1570 // generate AddNode()
1571 out << " " << GetPointerName() << "->AddNode";
1572 if (dnode->IsOverlapping()) out << "Overlap";
1573 out << "(" << dvol->GetPointerName() << ", " << icopy;
1574 if (!matrix->IsIdentity()) out << ", " << matrix->GetPointerName();
1575 out << ");" << std::endl;
1576 }
1577 // Recursive loop to daughters
1578 for (i=0; i<nd; i++) {
1579 dnode = GetNode(i);
1580 dvol = dnode->GetVolume();
1581 dvol->SavePrimitive(out,"d");
1582 }
1583 }
1584}
1585
1586////////////////////////////////////////////////////////////////////////////////
1587/// Reset SavePrimitive bits.
1588
1590{
1594}
1595
1596////////////////////////////////////////////////////////////////////////////////
1597/// Execute mouse actions on this volume.
1598
1600{
1602 if (!painter) return;
1603 painter->ExecuteVolumeEvent(this, event, px, py);
1604}
1605
1606////////////////////////////////////////////////////////////////////////////////
1607/// search a daughter inside the list of nodes
1608
1610{
1611 return ((TGeoNode*)fNodes->FindObject(name));
1612}
1613
1614////////////////////////////////////////////////////////////////////////////////
1615/// Get the index of a daughter within check_list by providing the node pointer.
1616
1617Int_t TGeoVolume::GetNodeIndex(const TGeoNode *node, Int_t *check_list, Int_t ncheck) const
1618{
1619 TGeoNode *current = 0;
1620 for (Int_t i=0; i<ncheck; i++) {
1621 current = (TGeoNode*)fNodes->At(check_list[i]);
1622 if (current==node) return check_list[i];
1623 }
1624 return -1;
1625}
1626
1627////////////////////////////////////////////////////////////////////////////////
1628/// get index number for a given daughter
1629
1631{
1632 TGeoNode *current = 0;
1633 Int_t nd = GetNdaughters();
1634 if (!nd) return -1;
1635 for (Int_t i=0; i<nd; i++) {
1636 current = (TGeoNode*)fNodes->At(i);
1637 if (current==node) return i;
1638 }
1639 return -1;
1640}
1641
1642////////////////////////////////////////////////////////////////////////////////
1643/// Get volume info for the browser.
1644
1646{
1647 TGeoVolume *vol = (TGeoVolume*)this;
1649 if (!painter) return 0;
1650 return (char*)painter->GetVolumeInfo(vol, px, py);
1651}
1652
1653////////////////////////////////////////////////////////////////////////////////
1654/// Returns true if cylindrical voxelization is optimal.
1655
1657{
1658 Int_t nd = GetNdaughters();
1659 if (!nd) return kFALSE;
1660 Int_t id;
1661 Int_t ncyl = 0;
1662 TGeoNode *node;
1663 for (id=0; id<nd; id++) {
1664 node = (TGeoNode*)fNodes->At(id);
1665 ncyl += node->GetOptimalVoxels();
1666 }
1667 if (ncyl>(nd/2)) return kTRUE;
1668 return kFALSE;
1669}
1670
1671////////////////////////////////////////////////////////////////////////////////
1672/// Provide a pointer name containing uid.
1673
1675{
1676 static TString name;
1677 name.Form("p%s_%zx", GetName(), (size_t)this);
1678 return name.Data();
1679}
1680
1681////////////////////////////////////////////////////////////////////////////////
1682/// Getter for optimization structure.
1683
1685{
1686 if (fVoxels && !fVoxels->IsInvalid()) return fVoxels;
1687 return NULL;
1688}
1689
1690////////////////////////////////////////////////////////////////////////////////
1691/// Move perspective view focus to this volume
1692
1694{
1696 if (painter) painter->GrabFocus();
1697}
1698
1699////////////////////////////////////////////////////////////////////////////////
1700/// Returns true if the volume is an assembly or a scaled assembly.
1701
1703{
1704 return fShape->IsAssembly();
1705}
1706
1707////////////////////////////////////////////////////////////////////////////////
1708/// Clone this volume.
1709/// build a volume with same name, shape and medium
1710
1712{
1713 TGeoVolume *vol = new TGeoVolume(GetName(), fShape, fMedium);
1714 Int_t i;
1715 // copy volume attributes
1716 vol->SetTitle(GetTitle());
1717 vol->SetLineColor(GetLineColor());
1718 vol->SetLineStyle(GetLineStyle());
1719 vol->SetLineWidth(GetLineWidth());
1720 vol->SetFillColor(GetFillColor());
1721 vol->SetFillStyle(GetFillStyle());
1722 // copy other attributes
1723 Int_t nbits = 8*sizeof(UInt_t);
1724 for (i=0; i<nbits; i++)
1725 vol->SetAttBit(1<<i, TGeoAtt::TestAttBit(1<<i));
1726 for (i=14; i<24; i++)
1727 vol->SetBit(1<<i, TestBit(1<<i));
1728
1729 // copy field
1730 vol->SetField(fField);
1731 // Set bits
1732 for (i=0; i<nbits; i++)
1733 vol->SetBit(1<<i, TObject::TestBit(1<<i));
1734 vol->SetBit(kVolumeClone);
1735 // copy nodes
1736// CloneNodesAndConnect(vol);
1737 vol->MakeCopyNodes(this);
1738 // if volume is divided, copy finder
1739 vol->SetFinder(fFinder);
1740 // copy voxels
1741 TGeoVoxelFinder *voxels = 0;
1742 if (fVoxels) {
1743 voxels = new TGeoVoxelFinder(vol);
1744 vol->SetVoxelFinder(voxels);
1745 }
1746 // copy option, uid
1747 vol->SetOption(fOption);
1748 vol->SetNumber(fNumber);
1749 vol->SetNtotal(fNtotal);
1750 // copy extensions
1754 return vol;
1755}
1756
1757////////////////////////////////////////////////////////////////////////////////
1758/// Clone the array of nodes.
1759
1761{
1762 if (!fNodes) return;
1763 TGeoNode *node;
1764 Int_t nd = fNodes->GetEntriesFast();
1765 if (!nd) return;
1766 // create new list of nodes
1767 TObjArray *list = new TObjArray(nd);
1768 // attach it to new volume
1769 newmother->SetNodes(list);
1770// ((TObject*)newmother)->SetBit(kVolumeImportNodes);
1771 for (Int_t i=0; i<nd; i++) {
1772 //create copies of nodes and add them to list
1773 node = GetNode(i)->MakeCopyNode();
1774 if (!node) {
1775 Fatal("CloneNodesAndConnect", "cannot make copy node");
1776 return;
1777 }
1778 node->SetMotherVolume(newmother);
1779 list->Add(node);
1780 }
1781}
1782
1783////////////////////////////////////////////////////////////////////////////////
1784/// make a new list of nodes and copy all nodes of other volume inside
1785
1787{
1788 Int_t nd = other->GetNdaughters();
1789 if (!nd) return;
1790 if (fNodes) {
1792 delete fNodes;
1793 }
1794 fNodes = new TObjArray();
1795 for (Int_t i=0; i<nd; i++) fNodes->Add(other->GetNode(i));
1797}
1798
1799////////////////////////////////////////////////////////////////////////////////
1800/// make a copy of this volume
1801/// build a volume with same name, shape and medium
1802
1804{
1805 TGeoVolume *vol = new TGeoVolume(GetName(), newshape, fMedium);
1806 // copy volume attributes
1807 vol->SetVisibility(IsVisible());
1808 vol->SetLineColor(GetLineColor());
1809 vol->SetLineStyle(GetLineStyle());
1810 vol->SetLineWidth(GetLineWidth());
1811 vol->SetFillColor(GetFillColor());
1812 vol->SetFillStyle(GetFillStyle());
1813 // copy field
1814 vol->SetField(fField);
1815 // if divided, copy division object
1816 if (fFinder) {
1817// Error("MakeCopyVolume", "volume %s divided", GetName());
1818 vol->SetFinder(fFinder);
1819 }
1820 // Copy extensions
1824// ((TObject*)vol)->SetBit(kVolumeImportNodes);
1825 ((TObject*)vol)->SetBit(kVolumeClone);
1827 return vol;
1828}
1829
1830////////////////////////////////////////////////////////////////////////////////
1831/// Make a copy of this volume which is reflected with respect to XY plane.
1832
1834{
1835 static TMap map(100);
1836 if (!fGeoManager->IsClosed()) {
1837 Error("MakeReflectedVolume", "Geometry must be closed.");
1838 return NULL;
1839 }
1840 TGeoVolume *vol = (TGeoVolume*)map.GetValue(this);
1841 if (vol) {
1842 if (newname && newname[0]) vol->SetName(newname);
1843 return vol;
1844 }
1845// printf("Making reflection for volume: %s\n", GetName());
1846 vol = CloneVolume();
1847 if (!vol) {
1848 Fatal("MakeReflectedVolume", "Cannot clone volume %s\n", GetName());
1849 return 0;
1850 }
1851 map.Add((TObject*)this, vol);
1852 if (newname && newname[0]) vol->SetName(newname);
1853 delete vol->GetNodes();
1854 vol->SetNodes(NULL);
1857 // The volume is now properly cloned, but with the same shape.
1858 // Reflect the shape (if any) and connect it.
1859 if (fShape) {
1860 TGeoShape *reflected_shape =
1862 vol->SetShape(reflected_shape);
1863 }
1864 // Reflect the daughters.
1865 Int_t nd = vol->GetNdaughters();
1866 if (!nd) return vol;
1867 TGeoNodeMatrix *node;
1868 TGeoMatrix *local, *local_cloned;
1869 TGeoVolume *new_vol;
1870 if (!vol->GetFinder()) {
1871 for (Int_t i=0; i<nd; i++) {
1872 node = (TGeoNodeMatrix*)vol->GetNode(i);
1873 local = node->GetMatrix();
1874// printf("%s before\n", node->GetName());
1875// local->Print();
1876 Bool_t reflected = local->IsReflection();
1877 local_cloned = new TGeoCombiTrans(*local);
1878 local_cloned->RegisterYourself();
1879 node->SetMatrix(local_cloned);
1880 if (!reflected) {
1881 // We need to reflect only the translation and propagate to daughters.
1882 // H' = Sz * H * Sz
1883 local_cloned->ReflectZ(kTRUE);
1884 local_cloned->ReflectZ(kFALSE);
1885// printf("%s after\n", node->GetName());
1886// node->GetMatrix()->Print();
1887 new_vol = node->GetVolume()->MakeReflectedVolume();
1888 node->SetVolume(new_vol);
1889 continue;
1890 }
1891 // The next daughter is already reflected, so reflect on Z everything and stop
1892 local_cloned->ReflectZ(kTRUE); // rot + tr
1893// printf("%s already reflected... After:\n", node->GetName());
1894// node->GetMatrix()->Print();
1895 }
1896 if (vol->GetVoxels()) vol->GetVoxels()->Voxelize();
1897 return vol;
1898 }
1899 // Volume is divided, so we have to reflect the division.
1900// printf(" ... divided %s\n", fFinder->ClassName());
1901 TGeoPatternFinder *new_finder = fFinder->MakeCopy(kTRUE);
1902 if (!new_finder) {
1903 Fatal("MakeReflectedVolume", "Could not copy finder for volume %s", GetName());
1904 return 0;
1905 }
1906 new_finder->SetVolume(vol);
1907 vol->SetFinder(new_finder);
1908 TGeoNodeOffset *nodeoff;
1909 new_vol = 0;
1910 for (Int_t i=0; i<nd; i++) {
1911 nodeoff = (TGeoNodeOffset*)vol->GetNode(i);
1912 nodeoff->SetFinder(new_finder);
1913 new_vol = nodeoff->GetVolume()->MakeReflectedVolume();
1914 nodeoff->SetVolume(new_vol);
1915 }
1916 return vol;
1917}
1918
1919////////////////////////////////////////////////////////////////////////////////
1920/// Set this volume as the TOP one (the whole geometry starts from here)
1921
1923{
1925}
1926
1927////////////////////////////////////////////////////////////////////////////////
1928/// Set the current tracking point.
1929
1931{
1933}
1934
1935////////////////////////////////////////////////////////////////////////////////
1936/// set the shape associated with this volume
1937
1939{
1940 if (!shape) {
1941 Error("SetShape", "No shape");
1942 return;
1943 }
1944 fShape = (TGeoShape*)shape;
1945}
1946
1947////////////////////////////////////////////////////////////////////////////////
1948/// sort nodes by decreasing volume of the bounding box. ONLY nodes comes first,
1949/// then overlapping nodes and finally division nodes.
1950
1952{
1953 if (!Valid()) {
1954 Error("SortNodes", "Bounding box not valid");
1955 return;
1956 }
1957 Int_t nd = GetNdaughters();
1958// printf("volume : %s, nd=%i\n", GetName(), nd);
1959 if (!nd) return;
1960 if (fFinder) return;
1961// printf("Nodes for %s\n", GetName());
1962 Int_t id = 0;
1963 TGeoNode *node = 0;
1964 TObjArray *nodes = new TObjArray(nd);
1965 Int_t inode = 0;
1966 // first put ONLY's
1967 for (id=0; id<nd; id++) {
1968 node = GetNode(id);
1969 if (node->InheritsFrom(TGeoNodeOffset::Class()) || node->IsOverlapping()) continue;
1970 nodes->Add(node);
1971// printf("inode %i ONLY\n", inode);
1972 inode++;
1973 }
1974 // second put overlapping nodes
1975 for (id=0; id<nd; id++) {
1976 node = GetNode(id);
1977 if (node->InheritsFrom(TGeoNodeOffset::Class()) || (!node->IsOverlapping())) continue;
1978 nodes->Add(node);
1979// printf("inode %i MANY\n", inode);
1980 inode++;
1981 }
1982 // third put the divided nodes
1983 if (fFinder) {
1984 fFinder->SetDivIndex(inode);
1985 for (id=0; id<nd; id++) {
1986 node = GetNode(id);
1987 if (!node->InheritsFrom(TGeoNodeOffset::Class())) continue;
1988 nodes->Add(node);
1989// printf("inode %i DIV\n", inode);
1990 inode++;
1991 }
1992 }
1993 if (inode != nd) printf(" volume %s : number of nodes does not match!!!\n", GetName());
1994 delete fNodes;
1995 fNodes = nodes;
1996}
1997
1998////////////////////////////////////////////////////////////////////////////////
1999/// Stream an object of class TGeoVolume.
2000
2002{
2003 if (R__b.IsReading()) {
2004 R__b.ReadClassBuffer(TGeoVolume::Class(), this);
2005 if (fVoxels && fVoxels->IsInvalid()) Voxelize("");
2006 } else {
2007 if (!fVoxels) {
2008 R__b.WriteClassBuffer(TGeoVolume::Class(), this);
2009 } else {
2011 TGeoVoxelFinder *voxels = fVoxels;
2012 fVoxels = 0;
2013 R__b.WriteClassBuffer(TGeoVolume::Class(), this);
2014 fVoxels = voxels;
2015 } else {
2016 R__b.WriteClassBuffer(TGeoVolume::Class(), this);
2017 }
2018 }
2019 }
2020}
2021
2022////////////////////////////////////////////////////////////////////////////////
2023/// Set the current options (none implemented)
2024
2026{
2027 fOption = option;
2028}
2029
2030////////////////////////////////////////////////////////////////////////////////
2031/// Set the line color.
2032
2034{
2035 TAttLine::SetLineColor(lcolor);
2036}
2037
2038////////////////////////////////////////////////////////////////////////////////
2039/// Set the line style.
2040
2042{
2043 TAttLine::SetLineStyle(lstyle);
2044}
2045
2046////////////////////////////////////////////////////////////////////////////////
2047/// Set the line width.
2048
2050{
2051 TAttLine::SetLineWidth(lwidth);
2052}
2053
2054////////////////////////////////////////////////////////////////////////////////
2055/// get the pointer to a daughter node
2056
2058{
2059 if (!fNodes) return 0;
2060 TGeoNode *node = (TGeoNode *)fNodes->FindObject(name);
2061 return node;
2062}
2063
2064////////////////////////////////////////////////////////////////////////////////
2065/// get the total size in bytes for this volume
2066
2068{
2069 Int_t count = 28+2+6+4+0; // TNamed+TGeoAtt+TAttLine+TAttFill+TAtt3D
2070 count += fName.Capacity() + fTitle.Capacity(); // name+title
2071 count += 7*sizeof(char*); // fShape + fMedium + fFinder + fField + fNodes + 2 extensions
2072 count += fOption.Capacity(); // fOption
2073 if (fShape) count += fShape->GetByteCount();
2074 if (fFinder) count += fFinder->GetByteCount();
2075 if (fNodes) {
2076 count += 32 + 4*fNodes->GetEntries(); // TObjArray
2077 TIter next(fNodes);
2078 TGeoNode *node;
2079 while ((node=(TGeoNode*)next())) count += node->GetByteCount();
2080 }
2081 return count;
2082}
2083
2084////////////////////////////////////////////////////////////////////////////////
2085/// loop all nodes marked as overlaps and find overlapping brothers
2086
2088{
2089 if (!Valid()) {
2090 Error("FindOverlaps","Bounding box not valid");
2091 return;
2092 }
2093 if (!fVoxels) return;
2094 Int_t nd = GetNdaughters();
2095 if (!nd) return;
2096 TGeoNode *node=0;
2097 Int_t inode = 0;
2098 for (inode=0; inode<nd; inode++) {
2099 node = GetNode(inode);
2100 if (!node->IsOverlapping()) continue;
2101 fVoxels->FindOverlaps(inode);
2102 }
2103}
2104
2105////////////////////////////////////////////////////////////////////////////////
2106/// Remove an existing daughter.
2107
2109{
2110 if (!fNodes || !fNodes->GetEntriesFast()) return;
2111 if (!fNodes->Remove(node)) return;
2112 fNodes->Compress();
2114 if (IsAssembly()) fShape->ComputeBBox();
2115}
2116
2117////////////////////////////////////////////////////////////////////////////////
2118/// Replace an existing daughter with a new volume having the same name but
2119/// possibly a new shape, position or medium. Not allowed for positioned assemblies.
2120/// For division cells, the new shape/matrix are ignored.
2121
2123{
2124 Int_t ind = GetIndex(nodeorig);
2125 if (ind < 0) return NULL;
2126 TGeoVolume *oldvol = nodeorig->GetVolume();
2127 if (oldvol->IsAssembly()) {
2128 Error("ReplaceNode", "Cannot replace node %s since it is an assembly", nodeorig->GetName());
2129 return NULL;
2130 }
2131 TGeoShape *shape = oldvol->GetShape();
2132 if (newshape && !nodeorig->IsOffset()) shape = newshape;
2133 TGeoMedium *med = oldvol->GetMedium();
2134 if (newmed) med = newmed;
2135 // Make a new volume
2136 TGeoVolume *vol = new TGeoVolume(oldvol->GetName(), shape, med);
2137 // copy volume attributes
2138 vol->SetVisibility(oldvol->IsVisible());
2139 vol->SetLineColor(oldvol->GetLineColor());
2140 vol->SetLineStyle(oldvol->GetLineStyle());
2141 vol->SetLineWidth(oldvol->GetLineWidth());
2142 vol->SetFillColor(oldvol->GetFillColor());
2143 vol->SetFillStyle(oldvol->GetFillStyle());
2144 // copy field
2145 vol->SetField(oldvol->GetField());
2146 // Make a copy of the node
2147 TGeoNode *newnode = nodeorig->MakeCopyNode();
2148 if (!newnode) {
2149 Fatal("ReplaceNode", "Cannot make copy node for %s", nodeorig->GetName());
2150 return 0;
2151 }
2152 // Change the volume for the new node
2153 newnode->SetVolume(vol);
2154 // Replace the matrix
2155 if (newpos && !nodeorig->IsOffset()) {
2156 TGeoNodeMatrix *nodemat = (TGeoNodeMatrix*)newnode;
2157 nodemat->SetMatrix(newpos);
2158 }
2159 // Replace nodeorig with new one
2160 fNodes->RemoveAt(ind);
2161 fNodes->AddAt(newnode, ind);
2163 if (IsAssembly()) fShape->ComputeBBox();
2164 return newnode;
2165}
2166
2167////////////////////////////////////////////////////////////////////////////////
2168/// Select this volume as matching an arbitrary criteria. The volume is added to
2169/// a static list and the flag TGeoVolume::kVolumeSelected is set. All flags need
2170/// to be reset at the end by calling the method with CLEAR=true. This will also clear
2171/// the list.
2172
2174{
2175 static TObjArray array(256);
2176 static Int_t len = 0;
2177 Int_t i;
2178 TObject *vol;
2179 if (clear) {
2180 for (i=0; i<len; i++) {
2181 vol = array.At(i);
2183 }
2184 array.Clear();
2185 len = 0;
2186 return;
2187 }
2189 array.AddAtAndExpand(this, len++);
2190}
2191
2192////////////////////////////////////////////////////////////////////////////////
2193/// set visibility of this volume
2194
2196{
2200 TSeqCollection *brlist = gROOT->GetListOfBrowsers();
2201 TIter next(brlist);
2202 TBrowser *browser = 0;
2203 while ((browser=(TBrowser*)next())) {
2204 browser->CheckObjectItem(this, vis);
2205 browser->Refresh();
2206 }
2207}
2208
2209////////////////////////////////////////////////////////////////////////////////
2210/// Set visibility for containers.
2211
2213{
2215 if (fGeoManager && fGeoManager->IsClosed()) {
2218 }
2219}
2220
2221////////////////////////////////////////////////////////////////////////////////
2222/// Set visibility for leaves.
2223
2225{
2227 if (fGeoManager && fGeoManager->IsClosed()) {
2230 }
2231}
2232
2233////////////////////////////////////////////////////////////////////////////////
2234/// Set visibility for leaves.
2235
2237{
2238 if (IsAssembly()) return;
2239 TGeoAtt::SetVisOnly(flag);
2240 if (fGeoManager && fGeoManager->IsClosed()) {
2243 }
2244}
2245
2246////////////////////////////////////////////////////////////////////////////////
2247/// Check if the shape of this volume is valid.
2248
2250{
2251 return fShape->IsValidBox();
2252}
2253
2254////////////////////////////////////////////////////////////////////////////////
2255/// Find a daughter node having VOL as volume and fill TGeoManager::fHMatrix
2256/// with its global matrix.
2257
2259{
2260 if (vol == this) return kTRUE;
2261 Int_t nd = GetNdaughters();
2262 if (!nd) return kFALSE;
2263 TGeoHMatrix *global = fGeoManager->GetHMatrix();
2264 if (!global) return kFALSE;
2265 TGeoNode *dnode;
2266 TGeoVolume *dvol;
2267 TGeoMatrix *local;
2268 Int_t i;
2269 for (i=0; i<nd; i++) {
2270 dnode = GetNode(i);
2271 dvol = dnode->GetVolume();
2272 if (dvol == vol) {
2273 local = dnode->GetMatrix();
2274 global->MultiplyLeft(local);
2275 return kTRUE;
2276 }
2277 }
2278 for (i=0; i<nd; i++) {
2279 dnode = GetNode(i);
2280 dvol = dnode->GetVolume();
2281 if (dvol->FindMatrixOfDaughterVolume(vol)) return kTRUE;
2282 }
2283 return kFALSE;
2284}
2285
2286////////////////////////////////////////////////////////////////////////////////
2287/// set visibility for daughters
2288
2290{
2291 SetVisDaughters(vis);
2294}
2295
2296////////////////////////////////////////////////////////////////////////////////
2297/// build the voxels for this volume
2298
2300{
2301 if (!Valid()) {
2302 Error("Voxelize", "Bounding box not valid");
2303 return;
2304 }
2305 // do not voxelize divided volumes
2306 if (fFinder) return;
2307 // or final leaves
2308 Int_t nd = GetNdaughters();
2309 if (!nd) return;
2310 // If this is an assembly, re-compute bounding box
2311 if (IsAssembly()) fShape->ComputeBBox();
2312 // delete old voxelization if any
2313 if (fVoxels) {
2315 fVoxels = 0;
2316 }
2317 // Create the voxels structure
2318 fVoxels = new TGeoVoxelFinder(this);
2320 if (fVoxels) {
2321 if (fVoxels->IsInvalid()) {
2322 delete fVoxels;
2323 fVoxels = 0;
2324 }
2325 }
2326}
2327
2328////////////////////////////////////////////////////////////////////////////////
2329/// Estimate the weight of a volume (in kg) with SIGMA(M)/M better than PRECISION.
2330/// Option can contain : v - verbose, a - analytical (default)
2331
2333{
2335 if (top != this) fGeoManager->SetTopVolume(this);
2336 else top = 0;
2337 Double_t weight = fGeoManager->Weight(precision, option);
2338 if (top) fGeoManager->SetTopVolume(top);
2339 return weight;
2340}
2341
2342////////////////////////////////////////////////////////////////////////////////
2343/// Analytical computation of the weight.
2344
2346{
2347 Double_t capacity = Capacity();
2348 Double_t weight = 0.0;
2349 Int_t i;
2350 Int_t nd = GetNdaughters();
2351 TGeoVolume *daughter;
2352 for (i=0; i<nd; i++) {
2353 daughter = GetNode(i)->GetVolume();
2354 weight += daughter->WeightA();
2355 capacity -= daughter->Capacity();
2356 }
2357 Double_t density = 0.0;
2358 if (!IsAssembly()) {
2359 if (fMedium) density = fMedium->GetMaterial()->GetDensity();
2360 if (density<0.01) density = 0.0; // do not weight gases
2361 }
2362 weight += 0.001*capacity * density; //[kg]
2363 return weight;
2364}
2365
2367
2368
2369////////////////////////////////////////////////////////////////////////////////
2370/// dummy constructor
2371
2373{
2374 fVolumes = 0;
2375 fDivision = 0;
2376 fNumed = 0;
2377 fNdiv = 0;
2378 fAxis = 0;
2379 fStart = 0;
2380 fStep = 0;
2381 fAttSet = kFALSE;
2383}
2384
2385////////////////////////////////////////////////////////////////////////////////
2386/// default constructor
2387
2389{
2390 fVolumes = new TObjArray();
2391 fDivision = 0;
2392 fNumed = 0;
2393 fNdiv = 0;
2394 fAxis = 0;
2395 fStart = 0;
2396 fStep = 0;
2397 fAttSet = kFALSE;
2399 SetName(name);
2400 SetMedium(med);
2401 fGeoManager->AddVolume(this);
2402// printf("--- volume multi %s created\n", name);
2403}
2404
2405////////////////////////////////////////////////////////////////////////////////
2406/// Destructor
2407
2409{
2410 if (fVolumes) delete fVolumes;
2411}
2412
2413////////////////////////////////////////////////////////////////////////////////
2414/// Add a volume with valid shape to the list of volumes. Copy all existing nodes
2415/// to this volume
2416
2418{
2419 Int_t idx = fVolumes->GetEntriesFast();
2420 fVolumes->AddAtAndExpand(vol,idx);
2421 vol->SetUniqueID(idx+1);
2422 TGeoVolumeMulti *div;
2423 TGeoVolume *cell;
2424 if (fDivision) {
2426 if (!div) {
2427 Fatal("AddVolume", "Cannot divide volume %s", vol->GetName());
2428 return;
2429 }
2430 for (Int_t i=0; i<div->GetNvolumes(); i++) {
2431 cell = div->GetVolume(i);
2432 fDivision->AddVolume(cell);
2433 }
2434 }
2435 if (fNodes) {
2436 Int_t nd = fNodes->GetEntriesFast();
2437 for (Int_t id=0; id<nd; id++) {
2438 TGeoNode *node = (TGeoNode*)fNodes->At(id);
2439 Bool_t many = node->IsOverlapping();
2440 if (many) vol->AddNodeOverlap(node->GetVolume(), node->GetNumber(), node->GetMatrix());
2441 else vol->AddNode(node->GetVolume(), node->GetNumber(), node->GetMatrix());
2442 }
2443 }
2444// vol->MakeCopyNodes(this);
2445}
2446
2447
2448////////////////////////////////////////////////////////////////////////////////
2449/// Add a new node to the list of nodes. This is the usual method for adding
2450/// daughters inside the container volume.
2451
2453{
2454 TGeoNode *n = TGeoVolume::AddNode(vol, copy_no, mat, option);
2455 Int_t nvolumes = fVolumes->GetEntriesFast();
2456 TGeoVolume *volume = 0;
2457 for (Int_t ivo=0; ivo<nvolumes; ivo++) {
2458 volume = GetVolume(ivo);
2459 volume->SetLineColor(GetLineColor());
2460 volume->SetLineStyle(GetLineStyle());
2461 volume->SetLineWidth(GetLineWidth());
2462 volume->SetVisibility(IsVisible());
2463 volume->AddNode(vol, copy_no, mat, option);
2464 }
2465 // printf("--- vmulti %s : node %s added to %i components\n", GetName(), vol->GetName(), nvolumes);
2466 return n;
2467}
2468
2469////////////////////////////////////////////////////////////////////////////////
2470/// Add a new node to the list of nodes, This node is possibly overlapping with other
2471/// daughters of the volume or extruding the volume.
2472
2474{
2475 TGeoVolume::AddNodeOverlap(vol, copy_no, mat, option);
2476 Int_t nvolumes = fVolumes->GetEntriesFast();
2477 TGeoVolume *volume = 0;
2478 for (Int_t ivo=0; ivo<nvolumes; ivo++) {
2479 volume = GetVolume(ivo);
2480 volume->SetLineColor(GetLineColor());
2481 volume->SetLineStyle(GetLineStyle());
2482 volume->SetLineWidth(GetLineWidth());
2483 volume->SetVisibility(IsVisible());
2484 volume->AddNodeOverlap(vol, copy_no, mat, option);
2485 }
2486// printf("--- vmulti %s : node ovlp %s added to %i components\n", GetName(), vol->GetName(), nvolumes);
2487}
2488
2489////////////////////////////////////////////////////////////////////////////////
2490/// Returns the last shape.
2491
2493{
2495 if (!vol) return 0;
2496 return vol->GetShape();
2497}
2498
2499////////////////////////////////////////////////////////////////////////////////
2500/// division of multiple volumes
2501
2502TGeoVolume *TGeoVolumeMulti::Divide(const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step, Int_t numed, const char *option)
2503{
2504 if (fDivision) {
2505 Error("Divide", "volume %s already divided", GetName());
2506 return 0;
2507 }
2508 Int_t nvolumes = fVolumes->GetEntriesFast();
2509 TGeoMedium *medium = fMedium;
2510 if (numed) {
2511 medium = fGeoManager->GetMedium(numed);
2512 if (!medium) {
2513 Error("Divide", "Invalid medium number %d for division volume %s", numed, divname);
2514 medium = fMedium;
2515 }
2516 }
2517 if (!nvolumes) {
2518 // this is a virtual volume
2519 fDivision = new TGeoVolumeMulti(divname, medium);
2520 fNumed = medium->GetId();
2521 fOption = option;
2522 fAxis = iaxis;
2523 fNdiv = ndiv;
2524 fStart = start;
2525 fStep = step;
2526 // nothing else to do at this stage
2527 return fDivision;
2528 }
2529 TGeoVolume *vol = 0;
2530 fDivision = new TGeoVolumeMulti(divname, medium);
2531 if (medium) fNumed = medium->GetId();
2532 fOption = option;
2533 fAxis = iaxis;
2534 fNdiv = ndiv;
2535 fStart = start;
2536 fStep = step;
2537 for (Int_t ivo=0; ivo<nvolumes; ivo++) {
2538 vol = GetVolume(ivo);
2539 vol->SetLineColor(GetLineColor());
2540 vol->SetLineStyle(GetLineStyle());
2541 vol->SetLineWidth(GetLineWidth());
2542 vol->SetVisibility(IsVisible());
2543 fDivision->AddVolume(vol->Divide(divname,iaxis,ndiv,start,step, numed, option));
2544 }
2545// printf("--- volume multi %s (%i volumes) divided\n", GetName(), nvolumes);
2546 if (numed) fDivision->SetMedium(medium);
2547 return fDivision;
2548}
2549
2550////////////////////////////////////////////////////////////////////////////////
2551/// Make a copy of this volume
2552/// build a volume with same name, shape and medium
2553
2555{
2556 TGeoVolume *vol = new TGeoVolume(GetName(), newshape, fMedium);
2557 Int_t i=0;
2558 // copy volume attributes
2559 vol->SetVisibility(IsVisible());
2560 vol->SetLineColor(GetLineColor());
2561 vol->SetLineStyle(GetLineStyle());
2562 vol->SetLineWidth(GetLineWidth());
2563 vol->SetFillColor(GetFillColor());
2564 vol->SetFillStyle(GetFillStyle());
2565 // copy field
2566 vol->SetField(fField);
2567 // Copy extensions
2570 // if divided, copy division object
2571// if (fFinder) {
2572// Error("MakeCopyVolume", "volume %s divided", GetName());
2573// vol->SetFinder(fFinder);
2574// }
2575 if (fDivision) {
2576 TGeoVolume *cell;
2578 if (!div) {
2579 Fatal("MakeCopyVolume", "Cannot divide volume %s", vol->GetName());
2580 return 0;
2581 }
2582 for (i=0; i<div->GetNvolumes(); i++) {
2583 cell = div->GetVolume(i);
2584 fDivision->AddVolume(cell);
2585 }
2586 }
2587
2588 if (!fNodes) return vol;
2589 TGeoNode *node;
2590 Int_t nd = fNodes->GetEntriesFast();
2591 if (!nd) return vol;
2592 // create new list of nodes
2593 TObjArray *list = new TObjArray();
2594 // attach it to new volume
2595 vol->SetNodes(list);
2596 ((TObject*)vol)->SetBit(kVolumeImportNodes);
2597 for (i=0; i<nd; i++) {
2598 //create copies of nodes and add them to list
2599 node = GetNode(i)->MakeCopyNode();
2600 if (!node) {
2601 Fatal("MakeCopyNode", "cannot make copy node for daughter %d of %s", i, GetName());
2602 return 0;
2603 }
2604 node->SetMotherVolume(vol);
2605 list->Add(node);
2606 }
2607 return vol;
2608}
2609
2610////////////////////////////////////////////////////////////////////////////////
2611/// Set the line color for all components.
2612
2614{
2616 Int_t nvolumes = fVolumes->GetEntriesFast();
2617 TGeoVolume *vol = 0;
2618 for (Int_t ivo=0; ivo<nvolumes; ivo++) {
2619 vol = GetVolume(ivo);
2620 vol->SetLineColor(lcolor);
2621 }
2622}
2623
2624////////////////////////////////////////////////////////////////////////////////
2625/// Set the line style for all components.
2626
2628{
2630 Int_t nvolumes = fVolumes->GetEntriesFast();
2631 TGeoVolume *vol = 0;
2632 for (Int_t ivo=0; ivo<nvolumes; ivo++) {
2633 vol = GetVolume(ivo);
2634 vol->SetLineStyle(lstyle);
2635 }
2636}
2637
2638////////////////////////////////////////////////////////////////////////////////
2639/// Set the line width for all components.
2640
2642{
2644 Int_t nvolumes = fVolumes->GetEntriesFast();
2645 TGeoVolume *vol = 0;
2646 for (Int_t ivo=0; ivo<nvolumes; ivo++) {
2647 vol = GetVolume(ivo);
2648 vol->SetLineWidth(lwidth);
2649 }
2650}
2651
2652////////////////////////////////////////////////////////////////////////////////
2653/// Set medium for a multiple volume.
2654
2656{
2658 Int_t nvolumes = fVolumes->GetEntriesFast();
2659 TGeoVolume *vol = 0;
2660 for (Int_t ivo=0; ivo<nvolumes; ivo++) {
2661 vol = GetVolume(ivo);
2662 vol->SetMedium(med);
2663 }
2664}
2665
2666
2667////////////////////////////////////////////////////////////////////////////////
2668/// Set visibility for all components.
2669
2671{
2673 Int_t nvolumes = fVolumes->GetEntriesFast();
2674 TGeoVolume *vol = 0;
2675 for (Int_t ivo=0; ivo<nvolumes; ivo++) {
2676 vol = GetVolume(ivo);
2677 vol->SetVisibility(vis);
2678 }
2679}
2680
2682
2683////////////////////////////////////////////////////////////////////////////////
2684/// Constructor.
2685
2687 fCurrent(-1), fNext(-1)
2688{
2689}
2690
2691////////////////////////////////////////////////////////////////////////////////
2692/// Destructor.
2693
2695{
2696}
2697
2698////////////////////////////////////////////////////////////////////////////////
2699
2701{
2703 return *fThreadData[tid];
2704}
2705
2706////////////////////////////////////////////////////////////////////////////////
2707
2709{
2710 std::lock_guard<std::mutex> guard(fMutex);
2712 std::vector<ThreadData_t*>::iterator i = fThreadData.begin();
2713 while (i != fThreadData.end())
2714 {
2715 delete *i;
2716 ++i;
2717 }
2718 fThreadData.clear();
2719 fThreadSize = 0;
2720}
2721
2722////////////////////////////////////////////////////////////////////////////////
2723
2725{
2726 std::lock_guard<std::mutex> guard(fMutex);
2727 // Create assembly thread data here
2728 fThreadData.resize(nthreads);
2729 fThreadSize = nthreads;
2730 for (Int_t tid=0; tid<nthreads; tid++) {
2731 if (fThreadData[tid] == 0) {
2732 fThreadData[tid] = new ThreadData_t;
2733 }
2734 }
2736}
2737
2738////////////////////////////////////////////////////////////////////////////////
2739
2741{
2742 return fThreadData[TGeoManager::ThreadId()]->fCurrent;
2743}
2744
2745////////////////////////////////////////////////////////////////////////////////
2746
2748{
2749 return fThreadData[TGeoManager::ThreadId()]->fNext;
2750}
2751
2752////////////////////////////////////////////////////////////////////////////////
2753
2755{
2757}
2758
2759////////////////////////////////////////////////////////////////////////////////
2760
2762{
2764}
2765
2766////////////////////////////////////////////////////////////////////////////////
2767/// Default constructor
2768
2770 :TGeoVolume()
2771{
2772 fThreadSize = 0;
2774}
2775
2776////////////////////////////////////////////////////////////////////////////////
2777/// Constructor. Just the name has to be provided. Assemblies does not have their own
2778/// shape or medium.
2779
2781 :TGeoVolume()
2782{
2783 fName = name;
2784 fName = fName.Strip();
2785 fShape = new TGeoShapeAssembly(this);
2787 fThreadSize = 0;
2789}
2790
2791////////////////////////////////////////////////////////////////////////////////
2792/// Destructor. The assembly is owner of its "shape".
2793
2795{
2797 if (fShape) delete fShape;
2798}
2799
2800////////////////////////////////////////////////////////////////////////////////
2801/// Add a component to the assembly.
2802
2804{
2805 TGeoNode *node = TGeoVolume::AddNode(vol, copy_no, mat, option);
2806 // ((TGeoShapeAssembly*)fShape)->RecomputeBoxLast();
2807 ((TGeoShapeAssembly*)fShape)->NeedsBBoxRecompute();
2808 return node;
2809}
2810
2811////////////////////////////////////////////////////////////////////////////////
2812/// Add an overlapping node - not allowed for assemblies.
2813
2815{
2816 Warning("AddNodeOverlap", "Declaring assembly %s as possibly overlapping inside %s not allowed. Using AddNode instead !",vol->GetName(),GetName());
2817 AddNode(vol, copy_no, mat, option);
2818}
2819
2820////////////////////////////////////////////////////////////////////////////////
2821/// Clone this volume.
2822/// build a volume with same name, shape and medium
2823
2825{
2827 Int_t i;
2828 // copy other attributes
2829 Int_t nbits = 8*sizeof(UInt_t);
2830 for (i=0; i<nbits; i++)
2831 vol->SetAttBit(1<<i, TGeoAtt::TestAttBit(1<<i));
2832 for (i=14; i<24; i++)
2833 vol->SetBit(1<<i, TestBit(1<<i));
2834
2835 // copy field
2836 vol->SetField(fField);
2837 // Set bits
2838 for (i=0; i<nbits; i++)
2839 vol->SetBit(1<<i, TObject::TestBit(1<<i));
2840 vol->SetBit(kVolumeClone);
2841 // make copy nodes
2842 vol->MakeCopyNodes(this);
2843// CloneNodesAndConnect(vol);
2844 ((TGeoShapeAssembly*)vol->GetShape())->NeedsBBoxRecompute();
2845 // copy voxels
2846 TGeoVoxelFinder *voxels = 0;
2847 if (fVoxels) {
2848 voxels = new TGeoVoxelFinder(vol);
2849 vol->SetVoxelFinder(voxels);
2850 }
2851 // copy option, uid
2852 vol->SetOption(fOption);
2853 vol->SetNumber(fNumber);
2854 vol->SetNtotal(fNtotal);
2855 vol->SetTitle(GetTitle());
2856 return vol;
2857}
2858
2859////////////////////////////////////////////////////////////////////////////////
2860/// Division makes no sense for assemblies.
2861
2863{
2864 Error("Divide","Assemblies cannot be divided");
2865 return 0;
2866}
2867
2868////////////////////////////////////////////////////////////////////////////////
2869/// Assign to the assembly a collection of identical volumes positioned according
2870/// a predefined pattern. The option can be spaced out or touching depending on the empty
2871/// space between volumes.
2872
2874{
2875 if (fNodes) {
2876 Error("Divide", "Cannot divide assembly %s since it has nodes", GetName());
2877 return NULL;
2878 }
2879 if (fFinder) {
2880 Error("Divide", "Assembly %s already divided", GetName());
2881 return NULL;
2882 }
2883 Int_t ncells = pattern->GetNdiv();
2884 if (!ncells || pattern->GetStep()<=0) {
2885 Error("Divide", "Pattern finder for dividing assembly %s not initialized. Use SetRange() method.", GetName());
2886 return NULL;
2887 }
2888 fFinder = pattern;
2889 TString opt(option);
2890 opt.ToLower();
2891 if (opt.Contains("spacedout")) fFinder->SetSpacedOut(kTRUE);
2893 // Position volumes
2894 for (Int_t i=0; i<ncells; i++) {
2895 fFinder->cd(i);
2896 TGeoNodeOffset *node = new TGeoNodeOffset(cell, i, 0.);
2897 node->SetFinder(fFinder);
2898 fNodes->Add(node);
2899 }
2900 return cell;
2901}
2902
2903////////////////////////////////////////////////////////////////////////////////
2904/// Make a clone of volume VOL but which is an assembly.
2905
2907{
2908 if (volorig->IsAssembly() || volorig->IsVolumeMulti()) return 0;
2909 Int_t nd = volorig->GetNdaughters();
2910 if (!nd) return 0;
2911 TGeoVolumeAssembly *vol = new TGeoVolumeAssembly(volorig->GetName());
2912 Int_t i;
2913 // copy other attributes
2914 Int_t nbits = 8*sizeof(UInt_t);
2915 for (i=0; i<nbits; i++)
2916 vol->SetAttBit(1<<i, volorig->TestAttBit(1<<i));
2917 for (i=14; i<24; i++)
2918 vol->SetBit(1<<i, volorig->TestBit(1<<i));
2919
2920 // copy field
2921 vol->SetField(volorig->GetField());
2922 // Set bits
2923 for (i=0; i<nbits; i++)
2924 vol->SetBit(1<<i, volorig->TestBit(1<<i));
2925 vol->SetBit(kVolumeClone);
2926 // make copy nodes
2927 vol->MakeCopyNodes(volorig);
2928// volorig->CloneNodesAndConnect(vol);
2929 vol->GetShape()->ComputeBBox();
2930 // copy voxels
2931 TGeoVoxelFinder *voxels = 0;
2932 if (volorig->GetVoxels()) {
2933 voxels = new TGeoVoxelFinder(vol);
2934 vol->SetVoxelFinder(voxels);
2935 }
2936 // copy option, uid
2937 vol->SetOption(volorig->GetOption());
2938 vol->SetNumber(volorig->GetNumber());
2939 vol->SetNtotal(volorig->GetNtotal());
2940 return vol;
2941}
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define s1(x)
Definition RSha256.hxx:91
short Style_t
Definition RtypesCore.h:89
int Int_t
Definition RtypesCore.h:45
short Color_t
Definition RtypesCore.h:92
unsigned int UInt_t
Definition RtypesCore.h:46
short Width_t
Definition RtypesCore.h:91
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
R__EXTERN TEnv * gEnv
Definition TEnv.h:170
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
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 char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t UChar_t len
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
R__EXTERN TGeoIdentity * gGeoIdentity
Definition TGeoMatrix.h:478
#define gROOT
Definition TROOT.h:405
R__EXTERN TStyle * gStyle
Definition TStyle.h:414
virtual Color_t GetFillColor() const
Return the fill area color.
Definition TAttFill.h:30
virtual Style_t GetFillStyle() const
Return the fill area style.
Definition TAttFill.h:31
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition TAttFill.h:39
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:33
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:42
virtual Width_t GetLineWidth() const
Return the line width.
Definition TAttLine.h:35
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
virtual Style_t GetLineStyle() const
Return the line style.
Definition TAttLine.h:34
virtual void SaveLineAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t widdef=1)
Save line attributes as C++ statement(s) on output stream out.
Definition TAttLine.cxx:273
Using a TBrowser one can browse all ROOT objects.
Definition TBrowser.h:37
void CheckObjectItem(TObject *obj, Bool_t check=kFALSE)
Change status of checkbox for this item.
Definition TBrowser.cxx:344
void Refresh()
Refresh browser contents.
Definition TBrowser.cxx:417
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=nullptr)=0
Bool_t IsReading() const
Definition TBuffer.h:86
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
TDirectory::TContext keeps track and restore the current directory.
Definition TDirectory.h:89
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
Definition TEnv.cxx:491
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition TFile.h:51
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:4053
Bool_t IsVisRaytrace() const
Definition TGeoAtt.h:87
virtual void SetVisOnly(Bool_t flag=kTRUE)
Set branch type visibility.
Definition TGeoAtt.cxx:95
Bool_t TestAttBit(UInt_t f) const
Definition TGeoAtt.h:68
virtual void SetVisLeaves(Bool_t flag=kTRUE)
Set branch type visibility.
Definition TGeoAtt.cxx:85
@ kSaveNodesAtt
Definition TGeoAtt.h:53
@ kSavePrimitiveAtt
Definition TGeoAtt.h:52
void SetVisDaughters(Bool_t vis=kTRUE)
Set visibility for the daughters.
Definition TGeoAtt.cxx:114
void ResetAttBit(UInt_t f)
Definition TGeoAtt.h:67
void SetVisRaytrace(Bool_t flag=kTRUE)
Definition TGeoAtt.h:70
Bool_t IsVisDaughters() const
Definition TGeoAtt.h:89
virtual void SetVisibility(Bool_t vis=kTRUE)
Set visibility for this object.
Definition TGeoAtt.cxx:105
void SetAttBit(UInt_t f)
Definition TGeoAtt.h:65
void SetVisTouched(Bool_t vis=kTRUE)
Mark visualization attributes as "modified".
Definition TGeoAtt.cxx:131
virtual void SetVisContainers(Bool_t flag=kTRUE)
Set branch type visibility.
Definition TGeoAtt.cxx:77
Class describing rotation + translation.
Definition TGeoMatrix.h:292
Composite shapes are Boolean combinations of two or more shape components.
void RegisterYourself()
Register the shape and all components to TGeoManager class.
ABC for user objects attached to TGeoVolume or TGeoNode.
virtual TGeoExtension * Grab()=0
virtual void Release() const =0
Matrix class used for computing global transformations Should NOT be used for node definition.
Definition TGeoMatrix.h:421
void MultiplyLeft(const TGeoMatrix *left)
multiply to the left with an other transformation if right is identity matrix, just return
The manager class for any TGeo geometry.
Definition TGeoManager.h:45
TObjArray * GetListOfOverlaps()
TList * GetListOfMedia() const
void SetUserPaintVolume(TGeoVolume *vol)
TObjArray * GetListOfVolumes() const
void ClearOverlaps()
Clear the list of overlaps.
TObjArray * GetListOfMatrices() const
Bool_t IsClosed() const
TVirtualGeoPainter * GetGeomPainter()
Make a default painter if none present. Returns pointer to it.
void SetVisOption(Int_t option=0)
set drawing mode :
Int_t AddMaterial(const TGeoMaterial *material)
Add a material to the list. Returns index of the material in list.
TGeoHMatrix * GetHMatrix()
Return stored current matrix (global matrix of the next touched node).
Int_t AddVolume(TGeoVolume *volume)
Add a volume to the list. Returns index of the volume in list.
Bool_t IsStreamingVoxels() const
void SetCurrentPoint(Double_t *point)
Int_t GetVisOption() const
Returns current depth to which geometry is drawn.
void SetTopVolume(TGeoVolume *vol)
Set the top volume and corresponding node as starting point of the geometry.
void ClearShape(const TGeoShape *shape)
Remove a shape from the list of shapes.
TGeoMedium * GetMedium(const char *medium) const
Search for a named tracking medium. All trailing blanks stripped.
Double_t Weight(Double_t precision=0.01, Option_t *option="va")
Estimate weight of volume VOL with a precision SIGMA(W)/W better than PRECISION.
void SetNsegments(Int_t nseg)
Set number of segments for approximating circles in drawing.
static UInt_t GetExportPrecision()
TVirtualGeoPainter * GetPainter() const
Bool_t IsCheckingOverlaps() const
void RandomPoints(const TGeoVolume *vol, Int_t npoints=10000, Option_t *option="")
Draw random points in the bounding box of a volume.
TList * GetListOfMaterials() const
void RandomRays(Int_t nrays=1000, Double_t startx=0, Double_t starty=0, Double_t startz=0, const char *target_vol=nullptr, Bool_t check_norm=kFALSE)
Randomly shoot nrays and plot intersections with surfaces for current top node.
void SetAllIndex()
Assigns uid's for all materials,media and matrices.
TObjArray * GetListOfShapes() const
TGeoVolume * GetTopVolume() const
static Int_t ThreadId()
Translates the current thread id to an ordinal number.
Int_t AddShape(const TGeoShape *shape)
Add a shape to the list. Returns index of the shape in list.
void SortOverlaps()
Sort overlaps by decreasing overlap distance. Extrusions comes first.
Base class describing materials.
virtual void Print(const Option_t *option="") const
print characteristics of this material
void SetUsed(Bool_t flag=kTRUE)
virtual Double_t GetDensity() const
Geometrical transformation package.
Definition TGeoMatrix.h:41
virtual void ReflectZ(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to XY.
Bool_t IsReflection() const
Definition TGeoMatrix.h:69
virtual void RegisterYourself()
Register the matrix in the current manager, which will become the owner.
void Print(Option_t *option="") const
print the matrix in 4x4 format
Bool_t IsIdentity() const
Definition TGeoMatrix.h:66
const char * GetPointerName() const
Provide a pointer name containing uid.
Bool_t IsRegistered() const
Definition TGeoMatrix.h:77
Media are used to store properties related to tracking and which are useful only when using geometry ...
Definition TGeoMedium.h:24
Int_t GetId() const
Definition TGeoMedium.h:48
TGeoMaterial * GetMaterial() const
Definition TGeoMedium.h:52
const char * GetPointerName() const
Provide a pointer name containing uid.
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
void SetMaterial(TGeoMaterial *mat)
Definition TGeoMedium.h:55
A node containing local transformation.
Definition TGeoNode.h:153
void SetMatrix(const TGeoMatrix *matrix)
Matrix setter.
Definition TGeoNode.cxx:767
virtual TGeoMatrix * GetMatrix() const
Definition TGeoNode.h:170
Node containing an offset.
Definition TGeoNode.h:184
void SetFinder(TGeoPatternFinder *finder)
Definition TGeoNode.h:206
static TClass * Class()
A node represent a volume positioned inside another.They store links to both volumes and to the TGeoM...
Definition TGeoNode.h:41
Bool_t IsOverlapping() const
Definition TGeoNode.h:105
TGeoVolume * GetVolume() const
Definition TGeoNode.h:97
void SetVolume(TGeoVolume *volume)
Definition TGeoNode.h:115
Bool_t IsOffset() const
Definition TGeoNode.h:103
virtual Int_t GetByteCount() const
Definition TGeoNode.h:84
void SetOverlapping(Bool_t flag=kTRUE)
Definition TGeoNode.h:118
virtual Int_t GetOptimalVoxels() const
Definition TGeoNode.h:99
virtual TGeoMatrix * GetMatrix() const =0
void SetMotherVolume(TGeoVolume *mother)
Definition TGeoNode.h:123
virtual TGeoNode * MakeCopyNode() const
Definition TGeoNode.h:111
void SetVirtual()
Definition TGeoNode.h:119
Int_t GetNumber() const
Definition TGeoNode.h:95
void SetNumber(Int_t number)
Definition TGeoNode.h:116
Base finder class for patterns.
virtual void cd(Int_t)
void SetSpacedOut(Bool_t flag)
void SetDivIndex(Int_t index)
virtual TGeoPatternFinder * MakeCopy(Bool_t reflect=kFALSE)=0
virtual Int_t GetDivAxis()
Int_t GetNdiv() const
void SetVolume(TGeoVolume *vol)
Double_t GetStep() const
void ClearThreadData() const
Double_t GetStart() const
virtual Int_t GetByteCount() const
void CreateThreadData(Int_t nthreads)
Create thread data for n threads max.
Class describing scale transformations.
Definition TGeoMatrix.h:245
static TGeoShape * MakeScaledShape(const char *name, TGeoShape *shape, TGeoScale *scale)
Create a scaled shape starting from a non-scaled one.
The shape encapsulating an assembly (union) of volumes.
Base abstract class for all shapes.
Definition TGeoShape.h:26
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const =0
virtual void CreateThreadData(Int_t)
Definition TGeoShape.h:68
Bool_t IsValid() const
Definition TGeoShape.h:140
virtual const char * GetAxisName(Int_t iaxis) const =0
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)=0
virtual Bool_t IsComposite() const
Definition TGeoShape.h:130
static Bool_t IsSameWithinTolerance(Double_t a, Double_t b)
Check if two numbers differ with less than a tolerance.
Bool_t IsRunTimeShape() const
Definition TGeoShape.h:139
virtual void ClearThreadData() const
Definition TGeoShape.h:67
const char * GetPointerName() const
Provide a pointer name containing uid.
void CheckShape(Int_t testNo, Int_t nsamples=10000, Option_t *option="")
Test for shape navigation methods.
virtual Bool_t IsValidBox() const =0
virtual const char * GetName() const
Get the shape name.
virtual Int_t GetByteCount() const =0
virtual void ComputeBBox()=0
virtual Double_t Capacity() const =0
@ kGeoSavePrimitive
Definition TGeoShape.h:65
virtual TGeoShape * GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix *mat) const =0
virtual Bool_t IsAssembly() const
Definition TGeoShape.h:129
Bool_t TestShapeBit(UInt_t f) const
Definition TGeoShape.h:163
Volume assemblies.
Definition TGeoVolume.h:305
static TGeoVolumeAssembly * MakeAssemblyFromVolume(TGeoVolume *vol)
Make a clone of volume VOL but which is an assembly.
virtual Int_t GetCurrentNodeIndex() const
virtual TGeoNode * AddNode(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat=nullptr, Option_t *option="")
Add a component to the assembly.
virtual TGeoVolume * Divide(const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step, Int_t numed=0, Option_t *option="")
Division makes no sense for assemblies.
virtual void ClearThreadData() const
virtual Int_t GetNextNodeIndex() const
virtual void AddNodeOverlap(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat, Option_t *option)
Add an overlapping node - not allowed for assemblies.
std::vector< ThreadData_t * > fThreadData
Definition TGeoVolume.h:321
virtual ~TGeoVolumeAssembly()
Destructor. The assembly is owner of its "shape".
std::mutex fMutex
Thread vector size.
Definition TGeoVolume.h:323
Int_t fThreadSize
Thread specific data vector.
Definition TGeoVolume.h:322
void SetNextNodeIndex(Int_t index)
virtual TGeoVolume * CloneVolume() const
Clone this volume.
void SetCurrentNodeIndex(Int_t index)
TGeoVolumeAssembly()
Default constructor.
virtual void CreateThreadData(Int_t nthreads)
ThreadData_t & GetThreadData() const
Volume families.
Definition TGeoVolume.h:255
virtual void SetVisibility(Bool_t vis=kTRUE)
Set visibility for all components.
virtual TGeoVolume * MakeCopyVolume(TGeoShape *newshape)
Make a copy of this volume build a volume with same name, shape and medium.
Double_t fStart
Definition TGeoVolume.h:262
void AddVolume(TGeoVolume *vol)
Add a volume with valid shape to the list of volumes.
virtual void SetLineStyle(Style_t lstyle)
Set the line style for all components.
virtual void AddNodeOverlap(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat, Option_t *option="")
Add a new node to the list of nodes, This node is possibly overlapping with other daughters of the vo...
TGeoVolume * GetVolume(Int_t id) const
Definition TGeoVolume.h:275
virtual TGeoNode * AddNode(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat, Option_t *option="")
Add a new node to the list of nodes.
virtual ~TGeoVolumeMulti()
Destructor.
virtual TGeoVolume * Divide(const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step, Int_t numed=0, Option_t *option="")
division of multiple volumes
virtual void SetLineColor(Color_t lcolor)
Set the line color for all components.
TGeoVolumeMulti * fDivision
Definition TGeoVolume.h:258
TGeoVolumeMulti()
dummy constructor
TGeoShape * GetLastShape() const
Returns the last shape.
Int_t GetNvolumes() const
Definition TGeoVolume.h:281
virtual void SetMedium(TGeoMedium *medium)
Set medium for a multiple volume.
TObjArray * fVolumes
Definition TGeoVolume.h:257
virtual void SetLineWidth(Width_t lwidth)
Set the line width for all components.
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:49
Double_t WeightA() const
Analytical computation of the weight.
void AddNodeOffset(TGeoVolume *vol, Int_t copy_no, Double_t offset=0, Option_t *option="")
Add a division node to the list of nodes.
virtual void cd(Int_t inode) const
Actualize matrix of node indexed <inode>
virtual void ClearThreadData() const
virtual ~TGeoVolume()
Destructor.
virtual void Print(Option_t *option="") const
Print volume info.
Bool_t IsVisContainers() const
Definition TGeoVolume.h:156
void SetVoxelFinder(TGeoVoxelFinder *finder)
Definition TGeoVolume.h:231
void RemoveNode(TGeoNode *node)
Remove an existing daughter.
Int_t GetNodeIndex(const TGeoNode *node, Int_t *check_list, Int_t ncheck) const
Get the index of a daughter within check_list by providing the node pointer.
Bool_t Valid() const
Check if the shape of this volume is valid.
Bool_t IsAllInvisible() const
Return TRUE if volume and all daughters are invisible.
Int_t fNtotal
Definition TGeoVolume.h:62
void MakeCopyNodes(const TGeoVolume *other)
make a new list of nodes and copy all nodes of other volume inside
void SetUserExtension(TGeoExtension *ext)
Connect user-defined extension to the volume.
TGeoExtension * GrabFWExtension() const
Get a copy of the framework extension pointer.
void SetNumber(Int_t number)
Definition TGeoVolume.h:233
void ClearNodes()
Definition TGeoVolume.h:100
void Raytrace(Bool_t flag=kTRUE)
Draw this volume with current settings and perform raytracing in the pad.
void RandomRays(Int_t nrays=10000, Double_t startx=0, Double_t starty=0, Double_t startz=0, const char *target_vol=nullptr, Bool_t check_norm=kFALSE)
Random raytracing method.
TGeoVolume()
dummy constructor
TGeoMedium * GetMedium() const
Definition TGeoVolume.h:174
virtual Bool_t IsFolder() const
Return TRUE if volume contains nodes.
void CloneNodesAndConnect(TGeoVolume *newmother) const
Clone the array of nodes.
void SortNodes()
sort nodes by decreasing volume of the bounding box.
Bool_t FindMatrixOfDaughterVolume(TGeoVolume *vol) const
Find a daughter node having VOL as volume and fill TGeoManager::fHMatrix with its global matrix.
void Voxelize(Option_t *option)
build the voxels for this volume
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
compute the closest distance of approach from point px,py to this volume
Double_t Capacity() const
Computes the capacity of this [cm^3] as the capacity of its shape.
virtual TGeoVolume * MakeCopyVolume(TGeoShape *newshape)
make a copy of this volume build a volume with same name, shape and medium
void ReplayCreation(const TGeoVolume *other)
Recreate the content of the other volume without pointer copying.
Double_t Weight(Double_t precision=0.01, Option_t *option="va")
Estimate the weight of a volume (in kg) with SIGMA(M)/M better than PRECISION.
Int_t fNumber
option - if any
Definition TGeoVolume.h:61
virtual void CreateThreadData(Int_t nthreads)
virtual Int_t GetByteCount() const
get the total size in bytes for this volume
virtual TGeoNode * AddNode(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat=nullptr, Option_t *option="")
Add a TGeoNode to the list of nodes.
Bool_t OptimizeVoxels()
Perform an extensive sampling to find which type of voxelization is most efficient.
void SetCurrentPoint(Double_t x, Double_t y, Double_t z)
Set the current tracking point.
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute mouse actions on this volume.
virtual void SetVisLeaves(Bool_t flag=kTRUE)
Set visibility for leaves.
void Browse(TBrowser *b)
How to browse a volume.
TGeoManager * fGeoManager
Definition TGeoVolume.h:57
TH2F * LegoPlot(Int_t ntheta=20, Double_t themin=0., Double_t themax=180., Int_t nphi=60, Double_t phimin=0., Double_t phimax=360., Double_t rmin=0., Double_t rmax=9999999, Option_t *option="")
Generate a lego plot fot the top volume, according to option.
TGeoVoxelFinder * fVoxels
Definition TGeoVolume.h:56
TGeoMaterial * GetMaterial() const
Definition TGeoVolume.h:173
virtual Bool_t IsVolumeMulti() const
Definition TGeoVolume.h:115
TGeoExtension * GrabUserExtension() const
Get a copy of the user extension pointer.
@ kVolumeSelected
Definition TGeoVolume.h:78
@ kVolumeImportNodes
Definition TGeoVolume.h:81
Int_t CountNodes(Int_t nlevels=1000, Int_t option=0)
Count total number of subnodes starting from this volume, nlevels down.
void GrabFocus()
Move perspective view focus to this volume.
void UnmarkSaved()
Reset SavePrimitive bits.
virtual TGeoVolume * CloneVolume() const
Clone this volume.
void SetFinder(TGeoPatternFinder *finder)
Definition TGeoVolume.h:232
Int_t GetNdaughters() const
Definition TGeoVolume.h:351
Bool_t IsValid() const
Definition TGeoVolume.h:153
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
void Grab()
Definition TGeoVolume.h:140
void CheckGeometry(Int_t nrays=1, Double_t startx=0, Double_t starty=0, Double_t startz=0) const
Shoot nrays with random directions from starting point (startx, starty, startz) in the reference fram...
void SelectVolume(Bool_t clear=kFALSE)
Select this volume as matching an arbitrary criteria.
const char * GetPointerName() const
Provide a pointer name containing uid.
static TClass * Class()
TObjArray * GetNodes()
Definition TGeoVolume.h:168
virtual char * GetObjectInfo(Int_t px, Int_t py) const
Get volume info for the browser.
Option_t * GetOption() const
Definition TGeoVolume.h:186
void ClearShape()
Clear the shape of this volume from the list held by the current manager.
void SetFWExtension(TGeoExtension *ext)
Connect framework defined extension to the volume.
void VisibleDaughters(Bool_t vis=kTRUE)
set visibility for daughters
void FindOverlaps() const
loop all nodes marked as overlaps and find overlapping brothers
virtual void AddNodeOverlap(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat=nullptr, Option_t *option="")
Add a TGeoNode to the list of nodes.
TGeoNode * GetNode(const char *name) const
get the pointer to a daughter node
virtual void SetVisibility(Bool_t vis=kTRUE)
set visibility of this volume
void RandomPoints(Int_t npoints=1000000, Option_t *option="")
Draw random points in the bounding box of this volume.
void CheckShapes()
check for negative parameters in shapes.
void SetNtotal(Int_t ntotal)
Definition TGeoVolume.h:234
virtual void Paint(Option_t *option="")
paint volume
Bool_t GetOptimalVoxels() const
Returns true if cylindrical voxelization is optimal.
TGeoNode * ReplaceNode(TGeoNode *nodeorig, TGeoShape *newshape=nullptr, TGeoMatrix *newpos=nullptr, TGeoMedium *newmed=nullptr)
Replace an existing daughter with a new volume having the same name but possibly a new shape,...
void InvisibleAll(Bool_t flag=kTRUE)
Make volume and each of it daughters (in)visible.
Bool_t IsVisibleDaughters() const
Definition TGeoVolume.h:155
TString fOption
just a hook for now
Definition TGeoVolume.h:60
void SaveAs(const char *filename, Option_t *option="") const
Save geometry having this as top volume as a C++ macro.
Int_t GetIndex(const TGeoNode *node) const
get index number for a given daughter
void SetNodes(TObjArray *nodes)
Definition TGeoVolume.h:215
TGeoPatternFinder * GetFinder() const
Definition TGeoVolume.h:176
void PrintVoxels() const
Print the voxels for this volume.
TGeoExtension * fUserExtension
Definition TGeoVolume.h:64
virtual void SetMedium(TGeoMedium *medium)
Definition TGeoVolume.h:230
TGeoVoxelFinder * GetVoxels() const
Getter for optimization structure.
void SetAttVisibility(Bool_t vis)
Definition TGeoVolume.h:221
void SetShape(const TGeoShape *shape)
set the shape associated with this volume
static TGeoMedium * DummyMedium()
TObject * fField
pointer to TGeoManager owning this volume
Definition TGeoVolume.h:59
Int_t GetNumber() const
Definition TGeoVolume.h:183
void CleanAll()
Clean data of the volume.
Bool_t IsTopVolume() const
True if this is the top volume of the geometry.
TGeoMedium * fMedium
Definition TGeoVolume.h:53
TGeoShape * GetShape() const
Definition TGeoVolume.h:189
void InspectMaterial() const
Inspect the material for this volume.
void PrintNodes() const
print nodes
virtual void Streamer(TBuffer &)
Stream an object of class TGeoVolume.
static TGeoMedium * fgDummyMedium
Definition TGeoVolume.h:54
void RegisterYourself(Option_t *option="")
Register the volume and all materials/media/matrices/shapes to the manager.
virtual void Draw(Option_t *option="")
draw top volume according to option
virtual void SetVisOnly(Bool_t flag=kTRUE)
Set visibility for leaves.
virtual TGeoVolume * Divide(const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step, Int_t numed=0, Option_t *option="")
Division a la G3.
Bool_t IsRaytracing() const
Check if the painter is currently ray-tracing the content of this volume.
TGeoShape * fShape
Definition TGeoVolume.h:52
void SetField(TObject *field)
Definition TGeoVolume.h:219
static TGeoVolume * Import(const char *filename, const char *name="", Option_t *option="")
Import a volume from a file.
Bool_t IsStyleDefault() const
check if the visibility and attributes are the default ones
static void CreateDummyMedium()
Create a dummy medium.
TGeoExtension * fFWExtension
Transient user-defined extension to volumes.
Definition TGeoVolume.h:65
void SetAsTopVolume()
Set this volume as the TOP one (the whole geometry starts from here)
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Bool_t IsVisLeaves() const
Definition TGeoVolume.h:157
virtual void SetVisContainers(Bool_t flag=kTRUE)
Set visibility for containers.
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
TObject * GetField() const
Definition TGeoVolume.h:175
TGeoPatternFinder * fFinder
dummy medium
Definition TGeoVolume.h:55
Int_t Export(const char *filename, const char *name="", Option_t *option="")
Export this volume to a file.
virtual void DrawOnly(Option_t *option="")
draw only this volume
virtual void SetLineColor(Color_t lcolor)
Set the line color.
void SetOption(const char *option)
Set the current options (none implemented)
virtual Bool_t IsAssembly() const
Returns true if the volume is an assembly or a scaled assembly.
TGeoVolume * MakeReflectedVolume(const char *newname="") const
Make a copy of this volume which is reflected with respect to XY plane.
TObjArray * fNodes
Definition TGeoVolume.h:51
virtual Bool_t IsVisible() const
Definition TGeoVolume.h:154
Int_t GetNtotal() const
Definition TGeoVolume.h:170
void InspectShape() const
Definition TGeoVolume.h:194
TGeoNode * FindNode(const char *name) const
search a daughter inside the list of nodes
void CheckOverlaps(Double_t ovlp=0.1, Option_t *option="") const
Overlap checking tool.
void SetOverlappingCandidate(Bool_t flag)
Definition TGeoVolume.h:216
Bool_t IsOverlappingCandidate() const
Definition TGeoVolume.h:147
Int_t fRefCount
Definition TGeoVolume.h:63
void CheckShape(Int_t testNo, Int_t nsamples=10000, Option_t *option="")
Tests for checking the shape navigation algorithms. See TGeoShape::CheckShape()
Finder class handling voxels.
Bool_t IsInvalid() const
void SetNeedRebuild(Bool_t flag=kTRUE)
virtual void Voxelize(Option_t *option="")
Voxelize attached volume according to option If the volume is an assembly, make sure the bbox is comp...
virtual void FindOverlaps(Int_t inode) const
create the list of nodes for which the bboxes overlap with inode's bbox
virtual void Print(Option_t *option="") const
Print the voxels.
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3060
2-D histogram with a float per channel (see TH1 documentation)}
Definition TH2.h:257
Book space in a file, create I/O buffers, to fill them, (un)compress them.
Definition TKey.h:28
virtual const char * GetClassName() const
Definition TKey.h:75
virtual TObject * ReadObj()
To read a TObject* from the file.
Definition TKey.cxx:750
TObject * FindObject(const char *name) const override
Find an object in this list using its name.
Definition TList.cxx:578
void Add(TObject *obj) override
Definition TList.h:81
TMap implements an associative array of (key,value) pairs using a THashTable for efficient retrieval ...
Definition TMap.h:40
void Add(TObject *obj) override
This function may not be used (but we need to provide it since it is a pure virtual in TCollection).
Definition TMap.cxx:54
TObject * GetValue(const char *keyname) const
Returns a pointer to the value associated with keyname as name of the key.
Definition TMap.cxx:236
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
TString fTitle
Definition TNamed.h:33
TString fName
Definition TNamed.h:32
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:140
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
void AddAt(TObject *obj, Int_t idx) override
Add object at position ids.
void Clear(Option_t *option="") override
Remove all objects from the array.
virtual void AddAtAndExpand(TObject *obj, Int_t idx)
Add object at position idx.
virtual void Compress()
Remove empty slots from array.
Int_t GetEntries() const override
Return the number of objects in array (i.e.
void Delete(Option_t *option="") override
Remove all objects from the array AND delete all heap based objects.
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
TObject * Remove(TObject *obj) override
Remove object from array.
TObject * RemoveAt(Int_t idx) override
Remove object at index idx.
TObject * FindObject(const char *name) const override
Find an object in this collection using its name.
void Add(TObject *obj) override
Definition TObjArray.h:68
Mother of all ROOT objects.
Definition TObject.h:41
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:201
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:207
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:956
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Definition TObject.cxx:745
virtual Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition TObject.cxx:874
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:774
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:525
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:970
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition TObject.cxx:998
virtual void SetUniqueID(UInt_t uid)
Set the unique object id.
Definition TObject.cxx:785
void ResetBit(UInt_t f)
Definition TObject.h:200
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:944
Sequenceable collection abstract base class.
Basic string class.
Definition TString.h:139
Ssiz_t Length() const
Definition TString.h:421
void ToLower()
Change string to lower-case.
Definition TString.cxx:1170
TSubString Strip(EStripType s=kTrailing, char c=' ') const
Return a substring of self stripped at beginning and/or end.
Definition TString.cxx:1151
const char * Data() const
Definition TString.h:380
Ssiz_t Capacity() const
Definition TString.h:368
Bool_t IsNull() const
Definition TString.h:418
TString & Remove(Ssiz_t pos)
Definition TString.h:685
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2356
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition TString.cxx:2334
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:636
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition TString.h:651
Abstract class for geometry painters.
virtual void SetTopVolume(TGeoVolume *vol)=0
virtual void ModifiedPad(Bool_t update=kFALSE) const =0
virtual TGeoVolume * GetTopVolume() const =0
virtual void DrawVolume(TGeoVolume *vol, Option_t *option="")=0
virtual void GrabFocus(Int_t nfr=0, Double_t dlong=0, Double_t dlat=0, Double_t dpsi=0)=0
virtual void CheckOverlaps(const TGeoVolume *vol, Double_t ovlp=0.1, Option_t *option="") const =0
virtual void Paint(Option_t *option="")=0
This method must be overridden if a class wants to paint itself.
virtual Int_t DistanceToPrimitiveVol(TGeoVolume *vol, Int_t px, Int_t py)=0
virtual Bool_t TestVoxels(TGeoVolume *vol)=0
virtual TGeoVolume * GetDrawnVolume() const =0
virtual const char * GetVolumeInfo(const TGeoVolume *volume, Int_t px, Int_t py) const =0
virtual void ExecuteVolumeEvent(TGeoVolume *volume, Int_t event, Int_t px, Int_t py)=0
virtual void CheckGeometry(Int_t nrays, Double_t startx, Double_t starty, Double_t startz) const =0
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
ThreadData_t()
index of next node to be entered