Logo ROOT  
Reference Guide
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 Geometry_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
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 "Riostream.h"
381#include "TString.h"
382#include "TBrowser.h"
383#include "TStyle.h"
384#include "TH2F.h"
385#include "TPad.h"
386#include "TROOT.h"
387#include "TClass.h"
388#include "TEnv.h"
389#include "TMap.h"
390#include "TFile.h"
391#include "TKey.h"
392
393#include "TGeoManager.h"
394#include "TGeoNode.h"
395#include "TGeoMatrix.h"
396#include "TVirtualGeoPainter.h"
397#include "TGeoVolume.h"
398#include "TGeoShapeAssembly.h"
399#include "TGeoScaledShape.h"
400#include "TGeoCompositeShape.h"
401#include "TGeoVoxelFinder.h"
402#include "TGeoExtension.h"
403
405
407
408////////////////////////////////////////////////////////////////////////////////
409/// Create a dummy medium
410
412{
413 if (fgDummyMedium) return;
415 fgDummyMedium->SetName("dummy");
416 TGeoMaterial *dummyMaterial = new TGeoMaterial();
417 dummyMaterial->SetName("dummy");
418 fgDummyMedium->SetMaterial(dummyMaterial);
419}
420
421////////////////////////////////////////////////////////////////////////////////
422
424{
427}
428
429////////////////////////////////////////////////////////////////////////////////
430
432{
433 if (fFinder) fFinder->CreateThreadData(nthreads);
434 if (fShape) fShape->CreateThreadData(nthreads);
435}
436
437////////////////////////////////////////////////////////////////////////////////
438
440{
441 return fgDummyMedium;
442}
443
444////////////////////////////////////////////////////////////////////////////////
445/// dummy constructor
446
448{
449 fNodes = 0;
450 fShape = 0;
451 fMedium = 0;
452 fFinder = 0;
453 fVoxels = 0;
455 fField = 0;
456 fOption = "";
457 fNumber = 0;
458 fNtotal = 0;
459 fRefCount = 0;
460 fUserExtension = 0;
461 fFWExtension = 0;
463}
464
465////////////////////////////////////////////////////////////////////////////////
466/// default constructor
467
468TGeoVolume::TGeoVolume(const char *name, const TGeoShape *shape, const TGeoMedium *med)
469 :TNamed(name, "")
470{
471 fName = fName.Strip();
472 fNodes = 0;
473 fShape = (TGeoShape*)shape;
474 if (fShape) {
476 Warning("Ctor", "volume %s has invalid shape", name);
477 }
478 if (!fShape->IsValid()) {
479 Fatal("ctor", "Shape of volume %s invalid. Aborting!", fName.Data());
480 }
481 }
482 fMedium = (TGeoMedium*)med;
484 fFinder = 0;
485 fVoxels = 0;
487 fField = 0;
488 fOption = "";
489 fNumber = 0;
490 fNtotal = 0;
491 fRefCount = 0;
492 fUserExtension = 0;
493 fFWExtension = 0;
496}
497
498////////////////////////////////////////////////////////////////////////////////
499/// Destructor
500
502{
503 if (fNodes) {
505 fNodes->Delete();
506 }
507 delete fNodes;
508 }
510 if (fVoxels) delete fVoxels;
513}
514
515////////////////////////////////////////////////////////////////////////////////
516/// How to browse a volume
517
519{
520 if (!b) return;
521
522// if (!GetNdaughters()) b->Add(this, GetName(), IsVisible());
523 TGeoVolume *daughter;
524 TString title;
525 for (Int_t i=0; i<GetNdaughters(); i++) {
526 daughter = GetNode(i)->GetVolume();
527 if(daughter->GetTitle()[0]) {
528 if (daughter->IsAssembly()) title.TString::Format("Assembly with %d daughter(s)",
529 daughter->GetNdaughters());
530 else if (daughter->GetFinder()) {
531 TString s1 = daughter->GetFinder()->ClassName();
532 s1.ReplaceAll("TGeoPattern","");
533 title.TString::Format("Volume having %s shape divided in %d %s slices",
534 daughter->GetShape()->ClassName(),daughter->GetNdaughters(), s1.Data());
535
536 } else title.TString::Format("Volume with %s shape having %d daughter(s)",
537 daughter->GetShape()->ClassName(),daughter->GetNdaughters());
538 daughter->SetTitle(title.Data());
539 }
540 b->Add(daughter, daughter->GetName(), daughter->IsVisible());
541// if (IsVisDaughters())
542// b->AddCheckBox(daughter, daughter->IsVisible());
543// else
544// b->AddCheckBox(daughter, kFALSE);
545 }
546}
547
548////////////////////////////////////////////////////////////////////////////////
549/// Computes the capacity of this [cm^3] as the capacity of its shape.
550/// In case of assemblies, the capacity is computed as the sum of daughter's capacities.
551
553{
554 if (!IsAssembly()) return fShape->Capacity();
555 Double_t capacity = 0.0;
556 Int_t nd = GetNdaughters();
557 Int_t i;
558 for (i=0; i<nd; i++) capacity += GetNode(i)->GetVolume()->Capacity();
559 return capacity;
560}
561
562////////////////////////////////////////////////////////////////////////////////
563/// Shoot nrays with random directions from starting point (startx, starty, startz)
564/// in the reference frame of this volume. Track each ray until exiting geometry, then
565/// shoot backwards from exiting point and compare boundary crossing points.
566
567void TGeoVolume::CheckGeometry(Int_t nrays, Double_t startx, Double_t starty, Double_t startz) const
568{
569 TGeoVolume *old_vol = fGeoManager->GetTopVolume();
570 if (old_vol!=this) fGeoManager->SetTopVolume((TGeoVolume*)this);
571 else old_vol=0;
574 painter->CheckGeometry(nrays, startx, starty, startz);
575}
576
577////////////////////////////////////////////////////////////////////////////////
578/// Overlap checking tool. Check for illegal overlaps within a limit OVLP.
579/// Use option="s[number]" to force overlap checking by sampling volume with
580/// [number] points.
581///
582/// Ex:
583/// ~~~ {.cpp}
584/// myVol->CheckOverlaps(0.01, "s10000000"); // shoot 10000000 points
585/// myVol->CheckOverlaps(0.01, "s"); // shoot the default value of 1e6 points
586/// ~~~
587
589{
590 if (!GetNdaughters() || fFinder) return;
591 Bool_t sampling = kFALSE;
592 TString opt(option);
593 opt.ToLower();
594 if (opt.Contains("s")) sampling = kTRUE;
596 if (!sampling) fGeoManager->SetNsegments(80);
599// Info("CheckOverlaps", "=== Checking overlaps for volume %s ===\n", GetName());
600 }
601 painter->CheckOverlaps(this, ovlp, option);
602// if (sampling) return;
606 Int_t novlps = overlaps->GetEntriesFast();
607 TNamed *obj;
609 for (Int_t i=0; i<novlps; i++) {
610 obj = (TNamed*)overlaps->At(i);
611 if (novlps<1000) name = TString::Format("ov%03d", i);
612 else name = TString::Format("ov%06d", i);
613 obj->SetName(name);
614 }
615 if (novlps) Info("CheckOverlaps", "Number of illegal overlaps/extrusions for volume %s: %d\n", GetName(), novlps);
616 }
617}
618
619////////////////////////////////////////////////////////////////////////////////
620/// Tests for checking the shape navigation algorithms. See TGeoShape::CheckShape()
621
622void TGeoVolume::CheckShape(Int_t testNo, Int_t nsamples, Option_t *option)
623{
624 fShape->CheckShape(testNo,nsamples,option);
625}
626
627////////////////////////////////////////////////////////////////////////////////
628/// Clean data of the volume.
629
631{
632 ClearNodes();
633 ClearShape();
634}
635
636////////////////////////////////////////////////////////////////////////////////
637/// Clear the shape of this volume from the list held by the current manager.
638
640{
642}
643
644////////////////////////////////////////////////////////////////////////////////
645/// check for negative parameters in shapes.
646
648{
649 if (fShape->IsRunTimeShape()) {
650 Error("CheckShapes", "volume %s has run-time shape", GetName());
651 InspectShape();
652 return;
653 }
654 if (!fNodes) return;
656 TGeoNode *node = 0;
657 TGeoNode *new_node;
658 const TGeoShape *shape = 0;
659 TGeoVolume *old_vol;
660 for (Int_t i=0; i<nd; i++) {
661 node=(TGeoNode*)fNodes->At(i);
662 // check if node has name
663 if (!node->GetName()[0]) printf("Daughter %i of volume %s - NO NAME!!!\n",
664 i, GetName());
665 old_vol = node->GetVolume();
666 shape = old_vol->GetShape();
667 if (shape->IsRunTimeShape()) {
668// printf(" Node %s/%s has shape with negative parameters. \n",
669// GetName(), node->GetName());
670// old_vol->InspectShape();
671 // make a copy of the node
672 new_node = node->MakeCopyNode();
673 if (!new_node) {
674 Fatal("CheckShapes", "Cannot make copy node for %s", node->GetName());
675 return;
676 }
677 TGeoShape *new_shape = shape->GetMakeRuntimeShape(fShape, node->GetMatrix());
678 if (!new_shape) {
679 Error("CheckShapes","cannot resolve runtime shape for volume %s/%s\n",
680 GetName(),old_vol->GetName());
681 continue;
682 }
683 TGeoVolume *new_volume = old_vol->MakeCopyVolume(new_shape);
684// printf(" new volume %s shape params :\n", new_volume->GetName());
685// new_volume->InspectShape();
686 new_node->SetVolume(new_volume);
687 // decouple the old node and put the new one instead
688 fNodes->AddAt(new_node, i);
689// new_volume->CheckShapes();
690 }
691 }
692}
693
694////////////////////////////////////////////////////////////////////////////////
695/// Count total number of subnodes starting from this volume, nlevels down
696/// - option = 0 (default) - count only once per volume
697/// - option = 1 - count every time
698/// - option = 2 - count volumes on visible branches
699/// - option = 3 - return maximum level counted already with option = 0
700
702{
703 static Int_t maxlevel = 0;
704 static Int_t nlev = 0;
705
706 if (option<0 || option>3) option = 0;
707 Int_t visopt = 0;
708 Int_t nd = GetNdaughters();
709 Bool_t last = (!nlevels || !nd)?kTRUE:kFALSE;
710 switch (option) {
711 case 0:
712 if (fNtotal) return fNtotal;
713 case 1:
714 fNtotal = 1;
715 break;
716 case 2:
717 visopt = fGeoManager->GetVisOption();
718 if (!IsVisDaughters()) last = kTRUE;
719 switch (visopt) {
721 fNtotal = (IsVisible())?1:0;
722 break;
724 fNtotal = (IsVisible() && last)?1:0;
725 }
726 if (!IsVisibleDaughters()) return fNtotal;
727 break;
728 case 3:
729 return maxlevel;
730 }
731 if (last) return fNtotal;
732 if (gGeoManager->GetTopVolume() == this) {
733 maxlevel=0;
734 nlev = 0;
735 }
736 if (nlev>maxlevel) maxlevel = nlev;
737 TGeoNode *node;
738 TGeoVolume *vol;
739 nlev++;
740 for (Int_t i=0; i<nd; i++) {
741 node = GetNode(i);
742 vol = node->GetVolume();
743 fNtotal += vol->CountNodes(nlevels-1, option);
744 }
745 nlev--;
746 return fNtotal;
747}
748
749////////////////////////////////////////////////////////////////////////////////
750/// Return TRUE if volume and all daughters are invisible.
751
753{
754 if (IsVisible()) return kFALSE;
755 Int_t nd = GetNdaughters();
756 for (Int_t i=0; i<nd; i++) if (GetNode(i)->GetVolume()->IsVisible()) return kFALSE;
757 return kTRUE;
758}
759
760////////////////////////////////////////////////////////////////////////////////
761/// Make volume and each of it daughters (in)visible.
762
764{
765 SetAttVisibility(!flag);
766 Int_t nd = GetNdaughters();
767 TObjArray *list = new TObjArray(nd+1);
768 list->Add(this);
769 TGeoVolume *vol;
770 for (Int_t i=0; i<nd; i++) {
771 vol = GetNode(i)->GetVolume();
772 vol->SetAttVisibility(!flag);
773 list->Add(vol);
774 }
775 TIter next(gROOT->GetListOfBrowsers());
776 TBrowser *browser = 0;
777 while ((browser=(TBrowser*)next())) {
778 for (Int_t i=0; i<nd+1; i++) {
779 vol = (TGeoVolume*)list->At(i);
780 browser->CheckObjectItem(vol, !flag);
781 }
782 browser->Refresh();
783 }
784 delete list;
786}
787
788////////////////////////////////////////////////////////////////////////////////
789/// Return TRUE if volume contains nodes
790
792{
793 return kTRUE;
794}
795
796////////////////////////////////////////////////////////////////////////////////
797/// check if the visibility and attributes are the default ones
798
800{
801 if (!IsVisible()) return kFALSE;
802 if (GetLineColor() != gStyle->GetLineColor()) return kFALSE;
803 if (GetLineStyle() != gStyle->GetLineStyle()) return kFALSE;
804 if (GetLineWidth() != gStyle->GetLineWidth()) return kFALSE;
805 return kTRUE;
806}
807
808////////////////////////////////////////////////////////////////////////////////
809/// True if this is the top volume of the geometry
810
812{
813 if (fGeoManager->GetTopVolume() == this) return kTRUE;
814 return kFALSE;
815}
816
817////////////////////////////////////////////////////////////////////////////////
818/// Check if the painter is currently ray-tracing the content of this volume.
819
821{
822 return TGeoAtt::IsVisRaytrace();
823}
824
825////////////////////////////////////////////////////////////////////////////////
826/// Inspect the material for this volume.
827
829{
830 GetMaterial()->Print();
831}
832
833////////////////////////////////////////////////////////////////////////////////
834/// Import a volume from a file.
835
836TGeoVolume *TGeoVolume::Import(const char *filename, const char *name, Option_t * /*option*/)
837{
838 if (!gGeoManager) gGeoManager = new TGeoManager("geometry","");
839 if (!filename) return 0;
840 TGeoVolume *volume = 0;
841 if (strstr(filename,".gdml")) {
842 // import from a gdml file
843 } else {
844 // import from a root file
846 TFile *f = TFile::Open(filename);
847 if (!f || f->IsZombie()) {
848 printf("Error: TGeoVolume::Import : Cannot open file %s\n", filename);
849 return 0;
850 }
851 if (name && name[0]) {
852 volume = (TGeoVolume*)f->Get(name);
853 } else {
854 TIter next(f->GetListOfKeys());
855 TKey *key;
856 while ((key = (TKey*)next())) {
857 if (strcmp(key->GetClassName(),"TGeoVolume") != 0) continue;
858 volume = (TGeoVolume*)key->ReadObj();
859 break;
860 }
861 }
862 delete f;
863 }
864 if (!volume) return NULL;
865 volume->RegisterYourself();
866 return volume;
867}
868
869////////////////////////////////////////////////////////////////////////////////
870/// Export this volume to a file.
871///
872/// - Case 1: root file or root/xml file
873/// if filename end with ".root". The key will be named name
874/// if filename end with ".xml" a root/xml file is produced.
875///
876/// - Case 2: C++ script
877/// if filename end with ".C"
878///
879/// - Case 3: gdml file
880/// if filename end with ".gdml"
881///
882/// NOTE that to use this option, the PYTHONPATH must be defined like
883/// export PYTHONPATH=$ROOTSYS/lib:$ROOTSYS/gdml
884///
885
886Int_t TGeoVolume::Export(const char *filename, const char *name, Option_t *option)
887{
888 TString sfile(filename);
889 if (sfile.Contains(".C")) {
890 //Save volume as a C++ script
891 Info("Export","Exporting volume %s as C++ code", GetName());
892 SaveAs(filename, "");
893 return 1;
894 }
895 if (sfile.Contains(".gdml")) {
896 //Save geometry as a gdml file
897 Info("Export","Exporting %s as gdml code - not implemented yet", GetName());
898 return 0;
899 }
900 if (sfile.Contains(".root") || sfile.Contains(".xml")) {
901 //Save volume in a root file
902 Info("Export","Exporting %s as root file.", GetName());
903 TString opt(option);
904 if (!opt.Length()) opt = "recreate";
905 TFile *f = TFile::Open(filename,opt.Data());
906 if (!f || f->IsZombie()) {
907 Error("Export","Cannot open file");
908 return 0;
909 }
910 TString keyname(name);
911 if (keyname.IsNull()) keyname = GetName();
912 Int_t nbytes = Write(keyname);
913 delete f;
914 return nbytes;
915 }
916 return 0;
917}
918
919////////////////////////////////////////////////////////////////////////////////
920/// Actualize matrix of node indexed <inode>
921
922void TGeoVolume::cd(Int_t inode) const
923{
924 if (fFinder) fFinder->cd(inode-fFinder->GetDivIndex());
925}
926
927////////////////////////////////////////////////////////////////////////////////
928/// Add a TGeoNode to the list of nodes. This is the usual method for adding
929/// daughters inside the container volume.
930
931void TGeoVolume::AddNode(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat, Option_t * /*option*/)
932{
933 TGeoMatrix *matrix = mat;
934 if (matrix==0) matrix = gGeoIdentity;
935 else matrix->RegisterYourself();
936 if (!vol) {
937 Error("AddNode", "Volume is NULL");
938 return;
939 }
940 if (!vol->IsValid()) {
941 Error("AddNode", "Won't add node with invalid shape");
942 printf("### invalid volume was : %s\n", vol->GetName());
943 return;
944 }
945 if (!fNodes) fNodes = new TObjArray();
946
947 if (fFinder) {
948 // volume already divided.
949 Error("AddNode", "Cannot add node %s_%i into divided volume %s", vol->GetName(), copy_no, GetName());
950 return;
951 }
952
953 TGeoNodeMatrix *node = 0;
954 node = new TGeoNodeMatrix(vol, matrix);
955 node->SetMotherVolume(this);
956 fNodes->Add(node);
957 TString name = TString::Format("%s_%d", vol->GetName(), copy_no);
958// if (fNodes->FindObject(name))
959// Warning("AddNode", "Volume %s : added node %s with same name", GetName(), name.Data());
960 node->SetName(name);
961 node->SetNumber(copy_no);
962 fRefCount++;
963 vol->Grab();
964}
965
966////////////////////////////////////////////////////////////////////////////////
967/// Add a division node to the list of nodes. The method is called by
968/// TGeoVolume::Divide() for creating the division nodes.
969
970void TGeoVolume::AddNodeOffset(TGeoVolume *vol, Int_t copy_no, Double_t offset, Option_t * /*option*/)
971{
972 if (!vol) {
973 Error("AddNodeOffset", "invalid volume");
974 return;
975 }
976 if (!vol->IsValid()) {
977 Error("AddNode", "Won't add node with invalid shape");
978 printf("### invalid volume was : %s\n", vol->GetName());
979 return;
980 }
981 if (!fNodes) fNodes = new TObjArray();
982 TGeoNode *node = new TGeoNodeOffset(vol, copy_no, offset);
983 node->SetMotherVolume(this);
984 fNodes->Add(node);
985 TString name = TString::Format("%s_%d", vol->GetName(), copy_no+1);
986 node->SetName(name);
987 node->SetNumber(copy_no+1);
988 vol->Grab();
989}
990
991////////////////////////////////////////////////////////////////////////////////
992/// Add a TGeoNode to the list of nodes. This is the usual method for adding
993/// daughters inside the container volume.
994
996{
997 if (!vol) {
998 Error("AddNodeOverlap", "Volume is NULL");
999 return;
1000 }
1001 if (!vol->IsValid()) {
1002 Error("AddNodeOverlap", "Won't add node with invalid shape");
1003 printf("### invalid volume was : %s\n", vol->GetName());
1004 return;
1005 }
1006 if (vol->IsAssembly()) {
1007 Warning("AddNodeOverlap", "Declaring assembly %s as possibly overlapping inside %s not allowed. Using AddNode instead !",vol->GetName(),GetName());
1008 AddNode(vol, copy_no, mat, option);
1009 return;
1010 }
1011 TGeoMatrix *matrix = mat;
1012 if (matrix==0) matrix = gGeoIdentity;
1013 else matrix->RegisterYourself();
1014 if (!fNodes) fNodes = new TObjArray();
1015
1016 if (fFinder) {
1017 // volume already divided.
1018 Error("AddNodeOverlap", "Cannot add node %s_%i into divided volume %s", vol->GetName(), copy_no, GetName());
1019 return;
1020 }
1021
1022 TGeoNodeMatrix *node = new TGeoNodeMatrix(vol, matrix);
1023 node->SetMotherVolume(this);
1024 fNodes->Add(node);
1025 TString name = TString::Format("%s_%d", vol->GetName(), copy_no);
1026 if (fNodes->FindObject(name))
1027 Warning("AddNode", "Volume %s : added node %s with same name", GetName(), name.Data());
1028 node->SetName(name);
1029 node->SetNumber(copy_no);
1030 node->SetOverlapping();
1031 if (vol->GetMedium() == fMedium)
1032 node->SetVirtual();
1033 vol->Grab();
1034}
1035
1036////////////////////////////////////////////////////////////////////////////////
1037/// Division a la G3. The volume will be divided along IAXIS (see shape classes), in NDIV
1038/// slices, from START with given STEP. The division volumes will have medium number NUMED.
1039/// If NUMED=0 they will get the medium number of the divided volume (this). If NDIV<=0,
1040/// all range of IAXIS will be divided and the resulting number of divisions will be centered on
1041/// IAXIS. If STEP<=0, the real STEP will be computed as the full range of IAXIS divided by NDIV.
1042/// Options (case insensitive):
1043/// - N - divide all range in NDIV cells (same effect as STEP<=0) (GSDVN in G3)
1044/// - NX - divide range starting with START in NDIV cells (GSDVN2 in G3)
1045/// - S - divide all range with given STEP. NDIV is computed and divisions will be centered
1046/// in full range (same effect as NDIV<=0) (GSDVS, GSDVT in G3)
1047/// - SX - same as DVS, but from START position. (GSDVS2, GSDVT2 in G3)
1048
1049TGeoVolume *TGeoVolume::Divide(const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step, Int_t numed, Option_t *option)
1050{
1051 if (fFinder) {
1052 // volume already divided.
1053 Fatal("Divide","volume %s already divided", GetName());
1054 return 0;
1055 }
1056 TString opt(option);
1057 opt.ToLower();
1058 TString stype = fShape->ClassName();
1059 if (!fNodes) fNodes = new TObjArray();
1060 Double_t xlo, xhi, range;
1061 range = fShape->GetAxisRange(iaxis, xlo, xhi);
1062 // for phi divisions correct the range
1063 if (!strcmp(fShape->GetAxisName(iaxis), "PHI")) {
1064 if ((start-xlo)<-1E-3) start+=360.;
1065 if (TGeoShape::IsSameWithinTolerance(range,360)) {
1066 xlo = start;
1067 xhi = start+range;
1068 }
1069 }
1070 if (range <=0) {
1071 InspectShape();
1072 Fatal("Divide", "cannot divide volume %s (%s) on %s axis", GetName(), stype.Data(), fShape->GetAxisName(iaxis));
1073 return 0;
1074 }
1075 if (ndiv<=0 || opt.Contains("s")) {
1076 if (step<=0) {
1077 Fatal("Divide", "invalid division type for volume %s : ndiv=%i, step=%g", GetName(), ndiv, step);
1078 return 0;
1079 }
1080 if (opt.Contains("x")) {
1081 if ((xlo-start)>1E-3 || (xhi-start)<-1E-3) {
1082 Fatal("Divide", "invalid START=%g for division on axis %s of volume %s. Range is (%g, %g)",
1083 start, fShape->GetAxisName(iaxis), GetName(), xlo, xhi);
1084 return 0;
1085 }
1086 xlo = start;
1087 range = xhi-xlo;
1088 }
1089 ndiv = Int_t((range+0.1*step)/step);
1090 Double_t ddx = range - ndiv*step;
1091 // always center the division in this case
1092 if (ddx>1E-3) Warning("Divide", "division of volume %s on %s axis (ndiv=%d) will be centered in the full range",
1093 GetName(), fShape->GetAxisName(iaxis), ndiv);
1094 start = xlo + 0.5*ddx;
1095 }
1096 if (step<=0 || opt.Contains("n")) {
1097 if (opt.Contains("x")) {
1098 if ((xlo-start)>1E-3 || (xhi-start)<-1E-3) {
1099 Fatal("Divide", "invalid START=%g for division on axis %s of volume %s. Range is (%g, %g)",
1100 start, fShape->GetAxisName(iaxis), GetName(), xlo, xhi);
1101 return 0;
1102 }
1103 xlo = start;
1104 range = xhi-xlo;
1105 }
1106 step = range/ndiv;
1107 start = xlo;
1108 }
1109
1110 Double_t end = start+ndiv*step;
1111 if (((start-xlo)<-1E-3) || ((end-xhi)>1E-3)) {
1112 Fatal("Divide", "division of volume %s on axis %s exceed range (%g, %g)",
1113 GetName(), fShape->GetAxisName(iaxis), xlo, xhi);
1114 return 0;
1115 }
1116 TGeoVolume *voldiv = fShape->Divide(this, divname, iaxis, ndiv, start, step);
1117 if (numed) {
1118 TGeoMedium *medium = fGeoManager->GetMedium(numed);
1119 if (!medium) {
1120 Fatal("Divide", "invalid medium number %d for division volume %s", numed, divname);
1121 return voldiv;
1122 }
1123 voldiv->SetMedium(medium);
1124 if (medium->GetMaterial()) medium->GetMaterial()->SetUsed();
1125 }
1126 return voldiv;
1127}
1128
1129////////////////////////////////////////////////////////////////////////////////
1130/// compute the closest distance of approach from point px,py to this volume
1131
1133{
1136 Int_t dist = 9999;
1137 if (!painter) return dist;
1138 dist = painter->DistanceToPrimitiveVol(this, px, py);
1139 return dist;
1140}
1141
1142////////////////////////////////////////////////////////////////////////////////
1143/// draw top volume according to option
1144
1146{
1151 if (!IsVisContainers()) SetVisLeaves();
1152 if (option && option[0] > 0) {
1153 painter->DrawVolume(this, option);
1154 } else {
1155 painter->DrawVolume(this, gEnv->GetValue("Viewer3D.DefaultDrawOption",""));
1156 }
1157}
1158
1159////////////////////////////////////////////////////////////////////////////////
1160/// draw only this volume
1161
1163{
1164 if (IsAssembly()) {
1165 Info("DrawOnly", "Volume assemblies do not support this option.");
1166 return;
1167 }
1169 SetVisOnly();
1172 if (option && option[0] > 0) {
1173 painter->DrawVolume(this, option);
1174 } else {
1175 painter->DrawVolume(this, gEnv->GetValue("Viewer3D.DefaultDrawOption",""));
1176 }
1177}
1178
1179////////////////////////////////////////////////////////////////////////////////
1180/// Perform an extensive sampling to find which type of voxelization is
1181/// most efficient.
1182
1184{
1185 printf("Optimizing volume %s ...\n", GetName());
1187 return painter->TestVoxels(this);
1188}
1189
1190////////////////////////////////////////////////////////////////////////////////
1191/// Print volume info
1192
1194{
1195 printf("== Volume: %s type %s positioned %d times\n", GetName(), ClassName(), fRefCount);
1196 InspectShape();
1198}
1199
1200////////////////////////////////////////////////////////////////////////////////
1201/// paint volume
1202
1204{
1206 painter->SetTopVolume(this);
1207// painter->Paint(option);
1208 if (option && option[0] > 0) {
1209 painter->Paint(option);
1210 } else {
1211 painter->Paint(gEnv->GetValue("Viewer3D.DefaultDrawOption",""));
1212 }
1213}
1214
1215////////////////////////////////////////////////////////////////////////////////
1216/// Print the voxels for this volume.
1217
1219{
1220 if (fVoxels) fVoxels->Print();
1221}
1222
1223////////////////////////////////////////////////////////////////////////////////
1224/// Recreate the content of the other volume without pointer copying. Voxels are
1225/// ignored and supposed to be created in a later step via Voxelize.
1226
1228{
1229 Int_t nd = other->GetNdaughters();
1230 if (!nd) return;
1231 TGeoPatternFinder *finder = other->GetFinder();
1232 if (finder) {
1233 Int_t iaxis = finder->GetDivAxis();
1234 Int_t ndiv = finder->GetNdiv();
1235 Double_t start = finder->GetStart();
1236 Double_t step = finder->GetStep();
1237 Int_t numed = other->GetNode(0)->GetVolume()->GetMedium()->GetId();
1238 TGeoVolume *voldiv = Divide(other->GetNode(0)->GetVolume()->GetName(), iaxis, ndiv, start, step, numed);
1239 voldiv->ReplayCreation(other->GetNode(0)->GetVolume());
1240 return;
1241 }
1242 for (Int_t i=0; i<nd; i++) {
1243 TGeoNode *node = other->GetNode(i);
1244 if (node->IsOverlapping()) AddNodeOverlap(node->GetVolume(), node->GetNumber(), node->GetMatrix());
1245 else AddNode(node->GetVolume(), node->GetNumber(), node->GetMatrix());
1246 }
1247}
1248
1249////////////////////////////////////////////////////////////////////////////////
1250/// print nodes
1251
1253{
1254 Int_t nd = GetNdaughters();
1255 for (Int_t i=0; i<nd; i++) {
1256 printf("%s\n", GetNode(i)->GetName());
1257 cd(i);
1258 GetNode(i)->GetMatrix()->Print();
1259 }
1260}
1261////////////////////////////////////////////////////////////////////////////////
1262/// Generate a lego plot fot the top volume, according to option.
1263
1265 Int_t nphi, Double_t phimin, Double_t phimax,
1266 Double_t rmin, Double_t rmax, Option_t *option)
1267{
1269 TGeoVolume *old_vol = fGeoManager->GetTopVolume();
1270 if (old_vol!=this) fGeoManager->SetTopVolume(this);
1271 else old_vol=0;
1272 TH2F *hist = p->LegoPlot(ntheta, themin, themax, nphi, phimin, phimax, rmin, rmax, option);
1273 hist->Draw("lego1sph");
1274 return hist;
1275}
1276
1277////////////////////////////////////////////////////////////////////////////////
1278/// Register the volume and all materials/media/matrices/shapes to the manager.
1279
1281{
1282 if (fGeoManager->GetListOfVolumes()->FindObject(this)) return;
1283 // Register volume
1284 fGeoManager->AddVolume(this);
1285 // Register shape
1287 if (fShape->IsComposite()) {
1289 comp->RegisterYourself();
1290 } else {
1292 }
1293 }
1294 // Register medium/material
1299 }
1300 // Register matrices for nodes.
1301 TGeoMatrix *matrix;
1302 TGeoNode *node;
1303 Int_t nd = GetNdaughters();
1304 Int_t i;
1305 for (i=0; i<nd; i++) {
1306 node = GetNode(i);
1307 matrix = node->GetMatrix();
1308 if (!matrix->IsRegistered()) matrix->RegisterYourself();
1309 else if (!fGeoManager->GetListOfMatrices()->FindObject(matrix)) {
1310 fGeoManager->GetListOfMatrices()->Add(matrix);
1311 }
1312 }
1313 // Call RegisterYourself recursively
1314 for (i=0; i<nd; i++) GetNode(i)->GetVolume()->RegisterYourself(option);
1315}
1316
1317////////////////////////////////////////////////////////////////////////////////
1318/// Draw random points in the bounding box of this volume.
1319
1321{
1323 TGeoVolume *old_vol = fGeoManager->GetTopVolume();
1324 if (old_vol!=this) fGeoManager->SetTopVolume(this);
1325 else old_vol=0;
1326 fGeoManager->RandomPoints(this, npoints, option);
1327 if (old_vol) fGeoManager->SetTopVolume(old_vol);
1328}
1329
1330////////////////////////////////////////////////////////////////////////////////
1331/// Random raytracing method.
1332
1333void TGeoVolume::RandomRays(Int_t nrays, Double_t startx, Double_t starty, Double_t startz, const char *target_vol, Bool_t check_norm)
1334{
1336 TGeoVolume *old_vol = fGeoManager->GetTopVolume();
1337 if (old_vol!=this) fGeoManager->SetTopVolume(this);
1338 else old_vol=0;
1339 fGeoManager->RandomRays(nrays, startx, starty, startz, target_vol, check_norm);
1340 if (old_vol) fGeoManager->SetTopVolume(old_vol);
1341}
1342
1343////////////////////////////////////////////////////////////////////////////////
1344/// Draw this volume with current settings and perform raytracing in the pad.
1345
1347{
1351 Bool_t drawn = (painter->GetDrawnVolume()==this)?kTRUE:kFALSE;
1352 if (!drawn) {
1353 painter->DrawVolume(this, "");
1355 painter->ModifiedPad();
1356 return;
1357 }
1359 painter->ModifiedPad();
1360}
1361
1362////////////////////////////////////////////////////////////////////////////////
1363/// Save geometry having this as top volume as a C++ macro.
1364
1365void TGeoVolume::SaveAs(const char *filename, Option_t *option) const
1366{
1367 if (!filename) return;
1368 std::ofstream out;
1369 out.open(filename, std::ios::out);
1370 if (out.bad()) {
1371 Error("SavePrimitive", "Bad file name: %s", filename);
1372 return;
1373 }
1375
1376 TString fname(filename);
1377 Int_t ind = fname.Index(".");
1378 if (ind>0) fname.Remove(ind);
1379 out << "void "<<fname<<"() {" << std::endl;
1380 out << " gSystem->Load(\"libGeom\");" << std::endl;
1382 out << std::setprecision(prec);
1383 ((TGeoVolume*)this)->SavePrimitive(out,option);
1384 out << "}" << std::endl;
1385}
1386
1387////////////////////////////////////////////////////////////////////////////////
1388/// Connect user-defined extension to the volume. The volume "grabs" a copy, so
1389/// the original object can be released by the producer. Release the previously
1390/// connected extension if any.
1391///
1392/// NOTE: This interface is intended for user extensions and is guaranteed not
1393/// to be used by TGeo
1394
1396{
1398 fUserExtension = 0;
1399 if (ext) fUserExtension = ext->Grab();
1400}
1401
1402////////////////////////////////////////////////////////////////////////////////
1403/// Connect framework defined extension to the volume. The volume "grabs" a copy,
1404/// so the original object can be released by the producer. Release the previously
1405/// connected extension if any.
1406///
1407/// NOTE: This interface is intended for the use by TGeo and the users should
1408/// NOT connect extensions using this method
1409
1411{
1413 fFWExtension = 0;
1414 if (ext) fFWExtension = ext->Grab();
1415}
1416
1417////////////////////////////////////////////////////////////////////////////////
1418/// Get a copy of the user extension pointer. The user must call Release() on
1419/// the copy pointer once this pointer is not needed anymore (equivalent to
1420/// delete() after calling new())
1421
1423{
1424 if (fUserExtension) return fUserExtension->Grab();
1425 return 0;
1426}
1427
1428////////////////////////////////////////////////////////////////////////////////
1429/// Get a copy of the framework extension pointer. The user must call Release() on
1430/// the copy pointer once this pointer is not needed anymore (equivalent to
1431/// delete() after calling new())
1432
1434{
1435 if (fFWExtension) return fFWExtension->Grab();
1436 return 0;
1437}
1438
1439////////////////////////////////////////////////////////////////////////////////
1440/// Save a primitive as a C++ statement(s) on output stream "out".
1441
1442void TGeoVolume::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
1443{
1444 Int_t i,icopy;
1445 Int_t nd = GetNdaughters();
1446 TGeoVolume *dvol;
1447 TGeoNode *dnode;
1448 TGeoMatrix *matrix;
1449
1450 // check if we need to save shape/volume
1451 Bool_t mustDraw = kFALSE;
1452 if (fGeoManager->GetGeomPainter()->GetTopVolume()==this) mustDraw = kTRUE;
1453 if (!option[0]) {
1455 out << " new TGeoManager(\"" << fGeoManager->GetName() << "\", \"" << fGeoManager->GetTitle() << "\");" << std::endl << std::endl;
1456// if (mustDraw) out << " Bool_t mustDraw = kTRUE;" << std::endl;
1457// else out << " Bool_t mustDraw = kFALSE;" << std::endl;
1458 out << " Double_t dx,dy,dz;" << std::endl;
1459 out << " Double_t dx1, dx2, dy1, dy2;" << std::endl;
1460 out << " Double_t vert[20], par[20];" << std::endl;
1461 out << " Double_t theta, phi, h1, bl1, tl1, alpha1, h2, bl2, tl2, alpha2;" << std::endl;
1462 out << " Double_t twist;" << std::endl;
1463 out << " Double_t origin[3];" << std::endl;
1464 out << " Double_t rmin, rmax, rmin1, rmax1, rmin2, rmax2;" << std::endl;
1465 out << " Double_t r, rlo, rhi;" << std::endl;
1466 out << " Double_t phi1, phi2;" << std::endl;
1467 out << " Double_t a,b;" << std::endl;
1468 out << " Double_t point[3], norm[3];" << std::endl;
1469 out << " Double_t rin, stin, rout, stout;" << std::endl;
1470 out << " Double_t thx, phx, thy, phy, thz, phz;" << std::endl;
1471 out << " Double_t alpha, theta1, theta2, phi1, phi2, dphi;" << std::endl;
1472 out << " Double_t tr[3], rot[9];" << std::endl;
1473 out << " Double_t z, density, radl, absl, w;" << std::endl;
1474 out << " Double_t lx,ly,lz,tx,ty,tz;" << std::endl;
1475 out << " Double_t xvert[50], yvert[50];" << std::endl;
1476 out << " Double_t zsect,x0,y0,scale0;" << std::endl;
1477 out << " Int_t nel, numed, nz, nedges, nvert;" << std::endl;
1478 out << " TGeoBoolNode *pBoolNode = 0;" << std::endl << std::endl;
1479 // first save materials/media
1480 out << " // MATERIALS, MIXTURES AND TRACKING MEDIA" << std::endl;
1481 SavePrimitive(out, "m");
1482 // then, save matrices
1483 out << std::endl << " // TRANSFORMATION MATRICES" << std::endl;
1484 SavePrimitive(out, "x");
1485 // save this volume and shape
1486 SavePrimitive(out, "s");
1487 out << std::endl << " // SET TOP VOLUME OF GEOMETRY" << std::endl;
1488 out << " gGeoManager->SetTopVolume(" << GetPointerName() << ");" << std::endl;
1489 // save daughters
1490 out << std::endl << " // SHAPES, VOLUMES AND GEOMETRICAL HIERARCHY" << std::endl;
1491 SavePrimitive(out, "d");
1492 out << std::endl << " // CLOSE GEOMETRY" << std::endl;
1493 out << " gGeoManager->CloseGeometry();" << std::endl;
1494 if (mustDraw) {
1495 if (!IsRaytracing()) out << " gGeoManager->GetTopVolume()->Draw();" << std::endl;
1496 else out << " gGeoManager->GetTopVolume()->Raytrace();" << std::endl;
1497 }
1498 return;
1499 }
1500 // check if we need to save shape/volume
1501 if (!strcmp(option, "s")) {
1502 // create the shape for this volume
1504 if (!IsAssembly()) {
1505 fShape->SavePrimitive(out,option);
1506 out << " // Volume: " << GetName() << std::endl;
1507 if (fMedium) out << " " << GetPointerName() << " = new TGeoVolume(\"" << GetName() << "\"," << fShape->GetPointerName() << ", "<< fMedium->GetPointerName() << ");" << std::endl;
1508 else out << " " << GetPointerName() << " = new TGeoVolume(\"" << GetName() << "\"," << fShape->GetPointerName() << ");" << std::endl;
1509
1510 } else {
1511 out << " // Assembly: " << GetName() << std::endl;
1512 out << " " << GetPointerName() << " = new TGeoVolumeAssembly(\"" << GetName() << "\"" << ");" << std::endl;
1513 }
1514 if (fLineColor != 1) out << " " << GetPointerName() << "->SetLineColor(" << fLineColor << ");" << std::endl;
1515 if (fLineWidth != 1) out << " " << GetPointerName() << "->SetLineWidth(" << fLineWidth << ");" << std::endl;
1516 if (fLineStyle != 1) out << " " << GetPointerName() << "->SetLineStyle(" << fLineStyle << ");" << std::endl;
1517 if (!IsVisible() && !IsAssembly()) out << " " << GetPointerName() << "->SetVisibility(kFALSE);" << std::endl;
1518 if (!IsVisibleDaughters()) out << " " << GetPointerName() << "->VisibleDaughters(kFALSE);" << std::endl;
1519 if (IsVisContainers()) out << " " << GetPointerName() << "->SetVisContainers(kTRUE);" << std::endl;
1520 if (IsVisLeaves()) out << " " << GetPointerName() << "->SetVisLeaves(kTRUE);" << std::endl;
1522 }
1523 // check if we need to save the media
1524 if (!strcmp(option, "m")) {
1525 if (fMedium) fMedium->SavePrimitive(out,option);
1526 for (i=0; i<nd; i++) {
1527 dvol = GetNode(i)->GetVolume();
1528 dvol->SavePrimitive(out,option);
1529 }
1530 return;
1531 }
1532 // check if we need to save the matrices
1533 if (!strcmp(option, "x")) {
1534 if (fFinder) {
1535 dvol = GetNode(0)->GetVolume();
1536 dvol->SavePrimitive(out,option);
1537 return;
1538 }
1539 for (i=0; i<nd; i++) {
1540 dnode = GetNode(i);
1541 matrix = dnode->GetMatrix();
1542 if (!matrix->IsIdentity()) matrix->SavePrimitive(out,option);
1543 dnode->GetVolume()->SavePrimitive(out,option);
1544 }
1545 return;
1546 }
1547 // check if we need to save volume daughters
1548 if (!strcmp(option, "d")) {
1549 if (!nd) return;
1552 if (fFinder) {
1553 // volume divided: generate volume->Divide()
1554 dnode = GetNode(0);
1555 dvol = dnode->GetVolume();
1556 out << " TGeoVolume *" << dvol->GetPointerName() << " = ";
1557 out << GetPointerName() << "->Divide(\"" << dvol->GetName() << "\", ";
1558 fFinder->SavePrimitive(out,option);
1559 if (fMedium != dvol->GetMedium()) {
1560 out << ", " << dvol->GetMedium()->GetId();
1561 }
1562 out << ");" << std::endl;
1563 dvol->SavePrimitive(out,"d");
1564 return;
1565 }
1566 for (i=0; i<nd; i++) {
1567 dnode = GetNode(i);
1568 dvol = dnode->GetVolume();
1569 dvol->SavePrimitive(out,"s");
1570 matrix = dnode->GetMatrix();
1571 icopy = dnode->GetNumber();
1572 // generate AddNode()
1573 out << " " << GetPointerName() << "->AddNode";
1574 if (dnode->IsOverlapping()) out << "Overlap";
1575 out << "(" << dvol->GetPointerName() << ", " << icopy;
1576 if (!matrix->IsIdentity()) out << ", " << matrix->GetPointerName();
1577 out << ");" << std::endl;
1578 }
1579 // Recursive loop to daughters
1580 for (i=0; i<nd; i++) {
1581 dnode = GetNode(i);
1582 dvol = dnode->GetVolume();
1583 dvol->SavePrimitive(out,"d");
1584 }
1585 }
1586}
1587
1588////////////////////////////////////////////////////////////////////////////////
1589/// Reset SavePrimitive bits.
1590
1592{
1596}
1597
1598////////////////////////////////////////////////////////////////////////////////
1599/// Execute mouse actions on this volume.
1600
1602{
1604 if (!painter) return;
1605 painter->ExecuteVolumeEvent(this, event, px, py);
1606}
1607
1608////////////////////////////////////////////////////////////////////////////////
1609/// search a daughter inside the list of nodes
1610
1612{
1613 return ((TGeoNode*)fNodes->FindObject(name));
1614}
1615
1616////////////////////////////////////////////////////////////////////////////////
1617/// Get the index of a daughter within check_list by providing the node pointer.
1618
1619Int_t TGeoVolume::GetNodeIndex(const TGeoNode *node, Int_t *check_list, Int_t ncheck) const
1620{
1621 TGeoNode *current = 0;
1622 for (Int_t i=0; i<ncheck; i++) {
1623 current = (TGeoNode*)fNodes->At(check_list[i]);
1624 if (current==node) return check_list[i];
1625 }
1626 return -1;
1627}
1628
1629////////////////////////////////////////////////////////////////////////////////
1630/// get index number for a given daughter
1631
1633{
1634 TGeoNode *current = 0;
1635 Int_t nd = GetNdaughters();
1636 if (!nd) return -1;
1637 for (Int_t i=0; i<nd; i++) {
1638 current = (TGeoNode*)fNodes->At(i);
1639 if (current==node) return i;
1640 }
1641 return -1;
1642}
1643
1644////////////////////////////////////////////////////////////////////////////////
1645/// Get volume info for the browser.
1646
1648{
1649 TGeoVolume *vol = (TGeoVolume*)this;
1651 if (!painter) return 0;
1652 return (char*)painter->GetVolumeInfo(vol, px, py);
1653}
1654
1655////////////////////////////////////////////////////////////////////////////////
1656/// Returns true if cylindrical voxelization is optimal.
1657
1659{
1660 Int_t nd = GetNdaughters();
1661 if (!nd) return kFALSE;
1662 Int_t id;
1663 Int_t ncyl = 0;
1664 TGeoNode *node;
1665 for (id=0; id<nd; id++) {
1666 node = (TGeoNode*)fNodes->At(id);
1667 ncyl += node->GetOptimalVoxels();
1668 }
1669 if (ncyl>(nd/2)) return kTRUE;
1670 return kFALSE;
1671}
1672
1673////////////////////////////////////////////////////////////////////////////////
1674/// Provide a pointer name containing uid.
1675
1677{
1678 static TString name;
1679 name = TString::Format("p%s_%lx", GetName(), (ULong_t)this);
1680 return (char*)name.Data();
1681}
1682
1683////////////////////////////////////////////////////////////////////////////////
1684/// Getter for optimization structure.
1685
1687{
1688 if (fVoxels && !fVoxels->IsInvalid()) return fVoxels;
1689 return NULL;
1690}
1691
1692////////////////////////////////////////////////////////////////////////////////
1693/// Move perspective view focus to this volume
1694
1696{
1698 if (painter) painter->GrabFocus();
1699}
1700
1701////////////////////////////////////////////////////////////////////////////////
1702/// Returns true if the volume is an assembly or a scaled assembly.
1703
1705{
1706 return fShape->IsAssembly();
1707}
1708
1709////////////////////////////////////////////////////////////////////////////////
1710/// Clone this volume.
1711/// build a volume with same name, shape and medium
1712
1714{
1715 TGeoVolume *vol = new TGeoVolume(GetName(), fShape, fMedium);
1716 Int_t i;
1717 // copy volume attributes
1718 vol->SetLineColor(GetLineColor());
1719 vol->SetLineStyle(GetLineStyle());
1720 vol->SetLineWidth(GetLineWidth());
1721 vol->SetFillColor(GetFillColor());
1722 vol->SetFillStyle(GetFillStyle());
1723 // copy other attributes
1724 Int_t nbits = 8*sizeof(UInt_t);
1725 for (i=0; i<nbits; i++)
1726 vol->SetAttBit(1<<i, TGeoAtt::TestAttBit(1<<i));
1727 for (i=14; i<24; i++)
1728 vol->SetBit(1<<i, TestBit(1<<i));
1729
1730 // copy field
1731 vol->SetField(fField);
1732 // Set bits
1733 for (i=0; i<nbits; i++)
1734 vol->SetBit(1<<i, TObject::TestBit(1<<i));
1735 vol->SetBit(kVolumeClone);
1736 // copy nodes
1737// CloneNodesAndConnect(vol);
1738 vol->MakeCopyNodes(this);
1739 // if volume is divided, copy finder
1740 vol->SetFinder(fFinder);
1741 // copy voxels
1742 TGeoVoxelFinder *voxels = 0;
1743 if (fVoxels) {
1744 voxels = new TGeoVoxelFinder(vol);
1745 vol->SetVoxelFinder(voxels);
1746 }
1747 // copy option, uid
1748 vol->SetOption(fOption);
1749 vol->SetNumber(fNumber);
1750 vol->SetNtotal(fNtotal);
1751 // copy extensions
1755 return vol;
1756}
1757
1758////////////////////////////////////////////////////////////////////////////////
1759/// Clone the array of nodes.
1760
1762{
1763 if (!fNodes) return;
1764 TGeoNode *node;
1765 Int_t nd = fNodes->GetEntriesFast();
1766 if (!nd) return;
1767 // create new list of nodes
1768 TObjArray *list = new TObjArray(nd);
1769 // attach it to new volume
1770 newmother->SetNodes(list);
1771// ((TObject*)newmother)->SetBit(kVolumeImportNodes);
1772 for (Int_t i=0; i<nd; i++) {
1773 //create copies of nodes and add them to list
1774 node = GetNode(i)->MakeCopyNode();
1775 if (!node) {
1776 Fatal("CloneNodesAndConnect", "cannot make copy node");
1777 return;
1778 }
1779 node->SetMotherVolume(newmother);
1780 list->Add(node);
1781 }
1782}
1783
1784////////////////////////////////////////////////////////////////////////////////
1785/// make a new list of nodes and copy all nodes of other volume inside
1786
1788{
1789 Int_t nd = other->GetNdaughters();
1790 if (!nd) return;
1791 if (fNodes) {
1793 delete fNodes;
1794 }
1795 fNodes = new TObjArray();
1796 for (Int_t i=0; i<nd; i++) fNodes->Add(other->GetNode(i));
1798}
1799
1800////////////////////////////////////////////////////////////////////////////////
1801/// make a copy of this volume
1802/// build a volume with same name, shape and medium
1803
1805{
1806 TGeoVolume *vol = new TGeoVolume(GetName(), newshape, fMedium);
1807 // copy volume attributes
1808 vol->SetVisibility(IsVisible());
1809 vol->SetLineColor(GetLineColor());
1810 vol->SetLineStyle(GetLineStyle());
1811 vol->SetLineWidth(GetLineWidth());
1812 vol->SetFillColor(GetFillColor());
1813 vol->SetFillStyle(GetFillStyle());
1814 // copy field
1815 vol->SetField(fField);
1816 // if divided, copy division object
1817 if (fFinder) {
1818// Error("MakeCopyVolume", "volume %s divided", GetName());
1819 vol->SetFinder(fFinder);
1820 }
1821 // Copy extensions
1825// ((TObject*)vol)->SetBit(kVolumeImportNodes);
1826 ((TObject*)vol)->SetBit(kVolumeClone);
1828 return vol;
1829}
1830
1831////////////////////////////////////////////////////////////////////////////////
1832/// Make a copy of this volume which is reflected with respect to XY plane.
1833
1835{
1836 static TMap map(100);
1837 if (!fGeoManager->IsClosed()) {
1838 Error("MakeReflectedVolume", "Geometry must be closed.");
1839 return NULL;
1840 }
1841 TGeoVolume *vol = (TGeoVolume*)map.GetValue(this);
1842 if (vol) {
1843 if (newname && newname[0]) vol->SetName(newname);
1844 return vol;
1845 }
1846// printf("Making reflection for volume: %s\n", GetName());
1847 vol = CloneVolume();
1848 if (!vol) {
1849 Fatal("MakeReflectedVolume", "Cannot clone volume %s\n", GetName());
1850 return 0;
1851 }
1852 map.Add((TObject*)this, vol);
1853 if (newname && newname[0]) vol->SetName(newname);
1854 delete vol->GetNodes();
1855 vol->SetNodes(NULL);
1858 // The volume is now properly cloned, but with the same shape.
1859 // Reflect the shape (if any) and connect it.
1860 if (fShape) {
1861 TGeoShape *reflected_shape =
1863 vol->SetShape(reflected_shape);
1864 }
1865 // Reflect the daughters.
1866 Int_t nd = vol->GetNdaughters();
1867 if (!nd) return vol;
1868 TGeoNodeMatrix *node;
1869 TGeoMatrix *local, *local_cloned;
1870 TGeoVolume *new_vol;
1871 if (!vol->GetFinder()) {
1872 for (Int_t i=0; i<nd; i++) {
1873 node = (TGeoNodeMatrix*)vol->GetNode(i);
1874 local = node->GetMatrix();
1875// printf("%s before\n", node->GetName());
1876// local->Print();
1877 Bool_t reflected = local->IsReflection();
1878 local_cloned = new TGeoCombiTrans(*local);
1879 local_cloned->RegisterYourself();
1880 node->SetMatrix(local_cloned);
1881 if (!reflected) {
1882 // We need to reflect only the translation and propagate to daughters.
1883 // H' = Sz * H * Sz
1884 local_cloned->ReflectZ(kTRUE);
1885 local_cloned->ReflectZ(kFALSE);
1886// printf("%s after\n", node->GetName());
1887// node->GetMatrix()->Print();
1888 new_vol = node->GetVolume()->MakeReflectedVolume();
1889 node->SetVolume(new_vol);
1890 continue;
1891 }
1892 // The next daughter is already reflected, so reflect on Z everything and stop
1893 local_cloned->ReflectZ(kTRUE); // rot + tr
1894// printf("%s already reflected... After:\n", node->GetName());
1895// node->GetMatrix()->Print();
1896 }
1897 if (vol->GetVoxels()) vol->GetVoxels()->Voxelize();
1898 return vol;
1899 }
1900 // Volume is divided, so we have to reflect the division.
1901// printf(" ... divided %s\n", fFinder->ClassName());
1902 TGeoPatternFinder *new_finder = fFinder->MakeCopy(kTRUE);
1903 if (!new_finder) {
1904 Fatal("MakeReflectedVolume", "Could not copy finder for volume %s", GetName());
1905 return 0;
1906 }
1907 new_finder->SetVolume(vol);
1908 vol->SetFinder(new_finder);
1909 TGeoNodeOffset *nodeoff;
1910 new_vol = 0;
1911 for (Int_t i=0; i<nd; i++) {
1912 nodeoff = (TGeoNodeOffset*)vol->GetNode(i);
1913 nodeoff->SetFinder(new_finder);
1914 new_vol = nodeoff->GetVolume()->MakeReflectedVolume();
1915 nodeoff->SetVolume(new_vol);
1916 }
1917 return vol;
1918}
1919
1920////////////////////////////////////////////////////////////////////////////////
1921/// Set this volume as the TOP one (the whole geometry starts from here)
1922
1924{
1926}
1927
1928////////////////////////////////////////////////////////////////////////////////
1929/// Set the current tracking point.
1930
1932{
1934}
1935
1936////////////////////////////////////////////////////////////////////////////////
1937/// set the shape associated with this volume
1938
1940{
1941 if (!shape) {
1942 Error("SetShape", "No shape");
1943 return;
1944 }
1945 fShape = (TGeoShape*)shape;
1946}
1947
1948////////////////////////////////////////////////////////////////////////////////
1949/// sort nodes by decreasing volume of the bounding box. ONLY nodes comes first,
1950/// then overlapping nodes and finally division nodes.
1951
1953{
1954 if (!Valid()) {
1955 Error("SortNodes", "Bounding box not valid");
1956 return;
1957 }
1958 Int_t nd = GetNdaughters();
1959// printf("volume : %s, nd=%i\n", GetName(), nd);
1960 if (!nd) return;
1961 if (fFinder) return;
1962// printf("Nodes for %s\n", GetName());
1963 Int_t id = 0;
1964 TGeoNode *node = 0;
1965 TObjArray *nodes = new TObjArray(nd);
1966 Int_t inode = 0;
1967 // first put ONLY's
1968 for (id=0; id<nd; id++) {
1969 node = GetNode(id);
1970 if (node->InheritsFrom(TGeoNodeOffset::Class()) || node->IsOverlapping()) continue;
1971 nodes->Add(node);
1972// printf("inode %i ONLY\n", inode);
1973 inode++;
1974 }
1975 // second put overlapping nodes
1976 for (id=0; id<nd; id++) {
1977 node = GetNode(id);
1978 if (node->InheritsFrom(TGeoNodeOffset::Class()) || (!node->IsOverlapping())) continue;
1979 nodes->Add(node);
1980// printf("inode %i MANY\n", inode);
1981 inode++;
1982 }
1983 // third put the divided nodes
1984 if (fFinder) {
1985 fFinder->SetDivIndex(inode);
1986 for (id=0; id<nd; id++) {
1987 node = GetNode(id);
1988 if (!node->InheritsFrom(TGeoNodeOffset::Class())) continue;
1989 nodes->Add(node);
1990// printf("inode %i DIV\n", inode);
1991 inode++;
1992 }
1993 }
1994 if (inode != nd) printf(" volume %s : number of nodes does not match!!!\n", GetName());
1995 delete fNodes;
1996 fNodes = nodes;
1997}
1998
1999////////////////////////////////////////////////////////////////////////////////
2000/// Stream an object of class TGeoVolume.
2001
2002void TGeoVolume::Streamer(TBuffer &R__b)
2003{
2004 if (R__b.IsReading()) {
2005 R__b.ReadClassBuffer(TGeoVolume::Class(), this);
2006 if (fVoxels && fVoxels->IsInvalid()) Voxelize("");
2007 } else {
2008 if (!fVoxels) {
2009 R__b.WriteClassBuffer(TGeoVolume::Class(), this);
2010 } else {
2012 TGeoVoxelFinder *voxels = fVoxels;
2013 fVoxels = 0;
2014 R__b.WriteClassBuffer(TGeoVolume::Class(), this);
2015 fVoxels = voxels;
2016 } else {
2017 R__b.WriteClassBuffer(TGeoVolume::Class(), this);
2018 }
2019 }
2020 }
2021}
2022
2023////////////////////////////////////////////////////////////////////////////////
2024/// Set the current options (none implemented)
2025
2026void TGeoVolume::SetOption(const char *option)
2027{
2028 fOption = option;
2029}
2030
2031////////////////////////////////////////////////////////////////////////////////
2032/// Set the line color.
2033
2035{
2036 TAttLine::SetLineColor(lcolor);
2037}
2038
2039////////////////////////////////////////////////////////////////////////////////
2040/// Set the line style.
2041
2043{
2044 TAttLine::SetLineStyle(lstyle);
2045}
2046
2047////////////////////////////////////////////////////////////////////////////////
2048/// Set the line width.
2049
2051{
2052 TAttLine::SetLineWidth(lwidth);
2053}
2054
2055////////////////////////////////////////////////////////////////////////////////
2056/// get the pointer to a daughter node
2057
2059{
2060 if (!fNodes) return 0;
2061 TGeoNode *node = (TGeoNode *)fNodes->FindObject(name);
2062 return node;
2063}
2064
2065////////////////////////////////////////////////////////////////////////////////
2066/// get the total size in bytes for this volume
2067
2069{
2070 Int_t count = 28+2+6+4+0; // TNamed+TGeoAtt+TAttLine+TAttFill+TAtt3D
2071 count += fName.Capacity() + fTitle.Capacity(); // name+title
2072 count += 7*sizeof(char*); // fShape + fMedium + fFinder + fField + fNodes + 2 extensions
2073 count += fOption.Capacity(); // fOption
2074 if (fShape) count += fShape->GetByteCount();
2075 if (fFinder) count += fFinder->GetByteCount();
2076 if (fNodes) {
2077 count += 32 + 4*fNodes->GetEntries(); // TObjArray
2078 TIter next(fNodes);
2079 TGeoNode *node;
2080 while ((node=(TGeoNode*)next())) count += node->GetByteCount();
2081 }
2082 return count;
2083}
2084
2085////////////////////////////////////////////////////////////////////////////////
2086/// loop all nodes marked as overlaps and find overlapping brothers
2087
2089{
2090 if (!Valid()) {
2091 Error("FindOverlaps","Bounding box not valid");
2092 return;
2093 }
2094 if (!fVoxels) return;
2095 Int_t nd = GetNdaughters();
2096 if (!nd) return;
2097 TGeoNode *node=0;
2098 Int_t inode = 0;
2099 for (inode=0; inode<nd; inode++) {
2100 node = GetNode(inode);
2101 if (!node->IsOverlapping()) continue;
2102 fVoxels->FindOverlaps(inode);
2103 }
2104}
2105
2106////////////////////////////////////////////////////////////////////////////////
2107/// Remove an existing daughter.
2108
2110{
2111 if (!fNodes || !fNodes->GetEntriesFast()) return;
2112 if (!fNodes->Remove(node)) return;
2113 fNodes->Compress();
2115 if (IsAssembly()) fShape->ComputeBBox();
2116}
2117
2118////////////////////////////////////////////////////////////////////////////////
2119/// Replace an existing daughter with a new volume having the same name but
2120/// possibly a new shape, position or medium. Not allowed for positioned assemblies.
2121/// For division cells, the new shape/matrix are ignored.
2122
2124{
2125 Int_t ind = GetIndex(nodeorig);
2126 if (ind < 0) return NULL;
2127 TGeoVolume *oldvol = nodeorig->GetVolume();
2128 if (oldvol->IsAssembly()) {
2129 Error("ReplaceNode", "Cannot replace node %s since it is an assembly", nodeorig->GetName());
2130 return NULL;
2131 }
2132 TGeoShape *shape = oldvol->GetShape();
2133 if (newshape && !nodeorig->IsOffset()) shape = newshape;
2134 TGeoMedium *med = oldvol->GetMedium();
2135 if (newmed) med = newmed;
2136 // Make a new volume
2137 TGeoVolume *vol = new TGeoVolume(oldvol->GetName(), shape, med);
2138 // copy volume attributes
2139 vol->SetVisibility(oldvol->IsVisible());
2140 vol->SetLineColor(oldvol->GetLineColor());
2141 vol->SetLineStyle(oldvol->GetLineStyle());
2142 vol->SetLineWidth(oldvol->GetLineWidth());
2143 vol->SetFillColor(oldvol->GetFillColor());
2144 vol->SetFillStyle(oldvol->GetFillStyle());
2145 // copy field
2146 vol->SetField(oldvol->GetField());
2147 // Make a copy of the node
2148 TGeoNode *newnode = nodeorig->MakeCopyNode();
2149 if (!newnode) {
2150 Fatal("ReplaceNode", "Cannot make copy node for %s", nodeorig->GetName());
2151 return 0;
2152 }
2153 // Change the volume for the new node
2154 newnode->SetVolume(vol);
2155 // Replace the matrix
2156 if (newpos && !nodeorig->IsOffset()) {
2157 TGeoNodeMatrix *nodemat = (TGeoNodeMatrix*)newnode;
2158 nodemat->SetMatrix(newpos);
2159 }
2160 // Replace nodeorig with new one
2161 fNodes->RemoveAt(ind);
2162 fNodes->AddAt(newnode, ind);
2164 if (IsAssembly()) fShape->ComputeBBox();
2165 return newnode;
2166}
2167
2168////////////////////////////////////////////////////////////////////////////////
2169/// Select this volume as matching an arbitrary criteria. The volume is added to
2170/// a static list and the flag TGeoVolume::kVolumeSelected is set. All flags need
2171/// to be reset at the end by calling the method with CLEAR=true. This will also clear
2172/// the list.
2173
2175{
2176 static TObjArray array(256);
2177 static Int_t len = 0;
2178 Int_t i;
2179 TObject *vol;
2180 if (clear) {
2181 for (i=0; i<len; i++) {
2182 vol = array.At(i);
2184 }
2185 array.Clear();
2186 len = 0;
2187 return;
2188 }
2190 array.AddAtAndExpand(this, len++);
2191}
2192
2193////////////////////////////////////////////////////////////////////////////////
2194/// set visibility of this volume
2195
2197{
2201 TSeqCollection *brlist = gROOT->GetListOfBrowsers();
2202 TIter next(brlist);
2203 TBrowser *browser = 0;
2204 while ((browser=(TBrowser*)next())) {
2205 browser->CheckObjectItem(this, vis);
2206 browser->Refresh();
2207 }
2208}
2209
2210////////////////////////////////////////////////////////////////////////////////
2211/// Set visibility for containers.
2212
2214{
2216 if (fGeoManager && fGeoManager->IsClosed()) {
2219 }
2220}
2221
2222////////////////////////////////////////////////////////////////////////////////
2223/// Set visibility for leaves.
2224
2226{
2228 if (fGeoManager && fGeoManager->IsClosed()) {
2231 }
2232}
2233
2234////////////////////////////////////////////////////////////////////////////////
2235/// Set visibility for leaves.
2236
2238{
2239 if (IsAssembly()) return;
2240 TGeoAtt::SetVisOnly(flag);
2241 if (fGeoManager && fGeoManager->IsClosed()) {
2244 }
2245}
2246
2247////////////////////////////////////////////////////////////////////////////////
2248/// Check if the shape of this volume is valid.
2249
2251{
2252 return fShape->IsValidBox();
2253}
2254
2255////////////////////////////////////////////////////////////////////////////////
2256/// Find a daughter node having VOL as volume and fill TGeoManager::fHMatrix
2257/// with its global matrix.
2258
2260{
2261 if (vol == this) return kTRUE;
2262 Int_t nd = GetNdaughters();
2263 if (!nd) return kFALSE;
2264 TGeoHMatrix *global = fGeoManager->GetHMatrix();
2265 if (!global) return kFALSE;
2266 TGeoNode *dnode;
2267 TGeoVolume *dvol;
2268 TGeoMatrix *local;
2269 Int_t i;
2270 for (i=0; i<nd; i++) {
2271 dnode = GetNode(i);
2272 dvol = dnode->GetVolume();
2273 if (dvol == vol) {
2274 local = dnode->GetMatrix();
2275 global->MultiplyLeft(local);
2276 return kTRUE;
2277 }
2278 }
2279 for (i=0; i<nd; i++) {
2280 dnode = GetNode(i);
2281 dvol = dnode->GetVolume();
2282 if (dvol->FindMatrixOfDaughterVolume(vol)) return kTRUE;
2283 }
2284 return kFALSE;
2285}
2286
2287////////////////////////////////////////////////////////////////////////////////
2288/// set visibility for daughters
2289
2291{
2292 SetVisDaughters(vis);
2295}
2296
2297////////////////////////////////////////////////////////////////////////////////
2298/// build the voxels for this volume
2299
2301{
2302 if (!Valid()) {
2303 Error("Voxelize", "Bounding box not valid");
2304 return;
2305 }
2306 // do not voxelize divided volumes
2307 if (fFinder) return;
2308 // or final leaves
2309 Int_t nd = GetNdaughters();
2310 if (!nd) return;
2311 // If this is an assembly, re-compute bounding box
2312 if (IsAssembly()) fShape->ComputeBBox();
2313 // delete old voxelization if any
2314 if (fVoxels) {
2316 fVoxels = 0;
2317 }
2318 // Create the voxels structure
2319 fVoxels = new TGeoVoxelFinder(this);
2320 fVoxels->Voxelize(option);
2321 if (fVoxels) {
2322 if (fVoxels->IsInvalid()) {
2323 delete fVoxels;
2324 fVoxels = 0;
2325 }
2326 }
2327}
2328
2329////////////////////////////////////////////////////////////////////////////////
2330/// Estimate the weight of a volume (in kg) with SIGMA(M)/M better than PRECISION.
2331/// Option can contain : v - verbose, a - analytical (default)
2332
2334{
2336 if (top != this) fGeoManager->SetTopVolume(this);
2337 else top = 0;
2338 Double_t weight = fGeoManager->Weight(precision, option);
2339 if (top) fGeoManager->SetTopVolume(top);
2340 return weight;
2341}
2342
2343////////////////////////////////////////////////////////////////////////////////
2344/// Analytical computation of the weight.
2345
2347{
2348 Double_t capacity = Capacity();
2349 Double_t weight = 0.0;
2350 Int_t i;
2351 Int_t nd = GetNdaughters();
2352 TGeoVolume *daughter;
2353 for (i=0; i<nd; i++) {
2354 daughter = GetNode(i)->GetVolume();
2355 weight += daughter->WeightA();
2356 capacity -= daughter->Capacity();
2357 }
2358 Double_t density = 0.0;
2359 if (!IsAssembly()) {
2360 if (fMedium) density = fMedium->GetMaterial()->GetDensity();
2361 if (density<0.01) density = 0.0; // do not weight gases
2362 }
2363 weight += 0.001*capacity * density; //[kg]
2364 return weight;
2365}
2366
2368
2369
2370////////////////////////////////////////////////////////////////////////////////
2371/// dummy constructor
2372
2374{
2375 fVolumes = 0;
2376 fDivision = 0;
2377 fNumed = 0;
2378 fNdiv = 0;
2379 fAxis = 0;
2380 fStart = 0;
2381 fStep = 0;
2382 fAttSet = kFALSE;
2384}
2385
2386////////////////////////////////////////////////////////////////////////////////
2387/// default constructor
2388
2390{
2391 fVolumes = new TObjArray();
2392 fDivision = 0;
2393 fNumed = 0;
2394 fNdiv = 0;
2395 fAxis = 0;
2396 fStart = 0;
2397 fStep = 0;
2398 fAttSet = kFALSE;
2400 SetName(name);
2401 SetMedium(med);
2402 fGeoManager->AddVolume(this);
2403// printf("--- volume multi %s created\n", name);
2404}
2405
2406////////////////////////////////////////////////////////////////////////////////
2407/// Destructor
2408
2410{
2411 if (fVolumes) delete fVolumes;
2412}
2413
2414////////////////////////////////////////////////////////////////////////////////
2415/// Add a volume with valid shape to the list of volumes. Copy all existing nodes
2416/// to this volume
2417
2419{
2420 Int_t idx = fVolumes->GetEntriesFast();
2421 fVolumes->AddAtAndExpand(vol,idx);
2422 vol->SetUniqueID(idx+1);
2423 TGeoVolumeMulti *div;
2424 TGeoVolume *cell;
2425 if (fDivision) {
2427 if (!div) {
2428 Fatal("AddVolume", "Cannot divide volume %s", vol->GetName());
2429 return;
2430 }
2431 for (Int_t i=0; i<div->GetNvolumes(); i++) {
2432 cell = div->GetVolume(i);
2433 fDivision->AddVolume(cell);
2434 }
2435 }
2436 if (fNodes) {
2437 Int_t nd = fNodes->GetEntriesFast();
2438 for (Int_t id=0; id<nd; id++) {
2439 TGeoNode *node = (TGeoNode*)fNodes->At(id);
2440 Bool_t many = node->IsOverlapping();
2441 if (many) vol->AddNodeOverlap(node->GetVolume(), node->GetNumber(), node->GetMatrix());
2442 else vol->AddNode(node->GetVolume(), node->GetNumber(), node->GetMatrix());
2443 }
2444 }
2445// vol->MakeCopyNodes(this);
2446}
2447
2448
2449////////////////////////////////////////////////////////////////////////////////
2450/// Add a new node to the list of nodes. This is the usual method for adding
2451/// daughters inside the container volume.
2452
2454{
2455 TGeoVolume::AddNode(vol, copy_no, mat, option);
2456 Int_t nvolumes = fVolumes->GetEntriesFast();
2457 TGeoVolume *volume = 0;
2458 for (Int_t ivo=0; ivo<nvolumes; ivo++) {
2459 volume = GetVolume(ivo);
2460 volume->SetLineColor(GetLineColor());
2461 volume->SetLineStyle(GetLineStyle());
2462 volume->SetLineWidth(GetLineWidth());
2463 volume->SetVisibility(IsVisible());
2464 volume->AddNode(vol, copy_no, mat, option);
2465 }
2466// printf("--- vmulti %s : node %s added to %i components\n", GetName(), vol->GetName(), nvolumes);
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{
2756 fThreadData[TGeoManager::ThreadId()]->fCurrent = index;
2757}
2758
2759////////////////////////////////////////////////////////////////////////////////
2760
2762{
2763 fThreadData[TGeoManager::ThreadId()]->fNext = index;
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 TGeoVolume::AddNode(vol,copy_no,mat,option);
2806// ((TGeoShapeAssembly*)fShape)->RecomputeBoxLast();
2807 ((TGeoShapeAssembly*)fShape)->NeedsBBoxRecompute();
2808}
2809
2810////////////////////////////////////////////////////////////////////////////////
2811/// Add an overlapping node - not allowed for assemblies.
2812
2814{
2815 Warning("AddNodeOverlap", "Declaring assembly %s as possibly overlapping inside %s not allowed. Using AddNode instead !",vol->GetName(),GetName());
2816 AddNode(vol, copy_no, mat, option);
2817}
2818
2819////////////////////////////////////////////////////////////////////////////////
2820/// Clone this volume.
2821/// build a volume with same name, shape and medium
2822
2824{
2826 Int_t i;
2827 // copy other attributes
2828 Int_t nbits = 8*sizeof(UInt_t);
2829 for (i=0; i<nbits; i++)
2830 vol->SetAttBit(1<<i, TGeoAtt::TestAttBit(1<<i));
2831 for (i=14; i<24; i++)
2832 vol->SetBit(1<<i, TestBit(1<<i));
2833
2834 // copy field
2835 vol->SetField(fField);
2836 // Set bits
2837 for (i=0; i<nbits; i++)
2838 vol->SetBit(1<<i, TObject::TestBit(1<<i));
2839 vol->SetBit(kVolumeClone);
2840 // make copy nodes
2841 vol->MakeCopyNodes(this);
2842// CloneNodesAndConnect(vol);
2843 ((TGeoShapeAssembly*)vol->GetShape())->NeedsBBoxRecompute();
2844 // copy voxels
2845 TGeoVoxelFinder *voxels = 0;
2846 if (fVoxels) {
2847 voxels = new TGeoVoxelFinder(vol);
2848 vol->SetVoxelFinder(voxels);
2849 }
2850 // copy option, uid
2851 vol->SetOption(fOption);
2852 vol->SetNumber(fNumber);
2853 vol->SetNtotal(fNtotal);
2854 return vol;
2855}
2856
2857////////////////////////////////////////////////////////////////////////////////
2858/// Division makes no sense for assemblies.
2859
2861{
2862 Error("Divide","Assemblies cannot be divided");
2863 return 0;
2864}
2865
2866////////////////////////////////////////////////////////////////////////////////
2867/// Assign to the assembly a collection of identical volumes positioned according
2868/// a predefined pattern. The option can be spaced out or touching depending on the empty
2869/// space between volumes.
2870
2872{
2873 if (fNodes) {
2874 Error("Divide", "Cannot divide assembly %s since it has nodes", GetName());
2875 return NULL;
2876 }
2877 if (fFinder) {
2878 Error("Divide", "Assembly %s already divided", GetName());
2879 return NULL;
2880 }
2881 Int_t ncells = pattern->GetNdiv();
2882 if (!ncells || pattern->GetStep()<=0) {
2883 Error("Divide", "Pattern finder for dividing assembly %s not initialized. Use SetRange() method.", GetName());
2884 return NULL;
2885 }
2886 fFinder = pattern;
2887 TString opt(option);
2888 opt.ToLower();
2889 if (opt.Contains("spacedout")) fFinder->SetSpacedOut(kTRUE);
2891 // Position volumes
2892 for (Int_t i=0; i<ncells; i++) {
2893 fFinder->cd(i);
2894 TGeoNodeOffset *node = new TGeoNodeOffset(cell, i, 0.);
2895 node->SetFinder(fFinder);
2896 fNodes->Add(node);
2897 }
2898 return cell;
2899}
2900
2901////////////////////////////////////////////////////////////////////////////////
2902/// Make a clone of volume VOL but which is an assembly.
2903
2905{
2906 if (volorig->IsAssembly() || volorig->IsVolumeMulti()) return 0;
2907 Int_t nd = volorig->GetNdaughters();
2908 if (!nd) return 0;
2909 TGeoVolumeAssembly *vol = new TGeoVolumeAssembly(volorig->GetName());
2910 Int_t i;
2911 // copy other attributes
2912 Int_t nbits = 8*sizeof(UInt_t);
2913 for (i=0; i<nbits; i++)
2914 vol->SetAttBit(1<<i, volorig->TestAttBit(1<<i));
2915 for (i=14; i<24; i++)
2916 vol->SetBit(1<<i, volorig->TestBit(1<<i));
2917
2918 // copy field
2919 vol->SetField(volorig->GetField());
2920 // Set bits
2921 for (i=0; i<nbits; i++)
2922 vol->SetBit(1<<i, volorig->TestBit(1<<i));
2923 vol->SetBit(kVolumeClone);
2924 // make copy nodes
2925 vol->MakeCopyNodes(volorig);
2926// volorig->CloneNodesAndConnect(vol);
2927 vol->GetShape()->ComputeBBox();
2928 // copy voxels
2929 TGeoVoxelFinder *voxels = 0;
2930 if (volorig->GetVoxels()) {
2931 voxels = new TGeoVoxelFinder(vol);
2932 vol->SetVoxelFinder(voxels);
2933 }
2934 // copy option, uid
2935 vol->SetOption(volorig->GetOption());
2936 vol->SetNumber(volorig->GetNumber());
2937 vol->SetNtotal(volorig->GetNtotal());
2938 return vol;
2939}
void Class()
Definition: Class.C:29
#define b(i)
Definition: RSha256.hxx:100
#define f(i)
Definition: RSha256.hxx:104
#define s1(x)
Definition: RSha256.hxx:91
int Int_t
Definition: RtypesCore.h:41
unsigned int UInt_t
Definition: RtypesCore.h:42
const Bool_t kFALSE
Definition: RtypesCore.h:88
unsigned long ULong_t
Definition: RtypesCore.h:51
short Width_t
Definition: RtypesCore.h:78
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
short Color_t
Definition: RtypesCore.h:79
short Style_t
Definition: RtypesCore.h:76
const Bool_t kTRUE
Definition: RtypesCore.h:87
const char Option_t
Definition: RtypesCore.h:62
#define ClassImp(name)
Definition: Rtypes.h:365
R__EXTERN TEnv * gEnv
Definition: TEnv.h:171
XFontStruct * id
Definition: TGX11.cxx:108
char name[80]
Definition: TGX11.cxx:109
R__EXTERN TGeoManager * gGeoManager
Definition: TGeoManager.h:601
R__EXTERN TGeoIdentity * gGeoIdentity
Definition: TGeoMatrix.h:478
#define gROOT
Definition: TROOT.h:415
R__EXTERN TStyle * gStyle
Definition: TStyle.h:407
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
Width_t fLineWidth
Line width.
Definition: TAttLine.h:23
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
Style_t fLineStyle
Line style.
Definition: TAttLine.h:22
Color_t fLineColor
Line color.
Definition: TAttLine.h:21
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:304
void Refresh()
Refresh browser contents.
Definition: TBrowser.cxx:377
Buffer base class used for serializing objects.
Definition: TBuffer.h:42
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
Bool_t IsReading() const
Definition: TBuffer.h:85
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
Small helper to keep current directory context.
Definition: TDirectory.h:41
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:48
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:3923
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
Class handling Boolean composition of shapes.
void RegisterYourself()
Register the shape and all components to TGeoManager class.
ABC for user objects attached to TGeoVolume or TGeoNode.
Definition: TGeoExtension.h:20
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:43
TObjArray * GetListOfOverlaps()
Definition: TGeoManager.h:488
TList * GetListOfMedia() const
Definition: TGeoManager.h:491
void SetUserPaintVolume(TGeoVolume *vol)
Definition: TGeoManager.h:237
TObjArray * GetListOfVolumes() const
Definition: TGeoManager.h:492
void ClearOverlaps()
Clear the list of overlaps.
TObjArray * GetListOfMatrices() const
Definition: TGeoManager.h:489
Bool_t IsClosed() const
Definition: TGeoManager.h:304
void RandomRays(Int_t nrays=1000, Double_t startx=0, Double_t starty=0, Double_t startz=0, const char *target_vol=0, Bool_t check_norm=kFALSE)
Randomly shoot nrays and plot intersections with surfaces for current top node.
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
Definition: TGeoManager.h:482
void SetCurrentPoint(Double_t *point)
Definition: TGeoManager.h:534
Int_t GetVisOption() const
Returns current depth to which geometry is drawn.
static UInt_t GetExportPrecision()
Definition: TGeoManager.h:478
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.
TVirtualGeoPainter * GetPainter() const
Definition: TGeoManager.h:212
Bool_t IsCheckingOverlaps() const
Definition: TGeoManager.h:410
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
Definition: TGeoManager.h:490
void SetAllIndex()
Assigns uid's for all materials,media and matrices.
TObjArray * GetListOfShapes() const
Definition: TGeoManager.h:494
TGeoVolume * GetTopVolume() const
Definition: TGeoManager.h:531
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.
Definition: TGeoMaterial.h:31
virtual void Print(const Option_t *option="") const
print characteristics of this material
void SetUsed(Bool_t flag=kTRUE)
Definition: TGeoMaterial.h:133
virtual Double_t GetDensity() const
Definition: TGeoMaterial.h:103
Geometrical transformation package.
Definition: TGeoMatrix.h:41
virtual void ReflectZ(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to XY.
Definition: TGeoMatrix.cxx:519
Bool_t IsReflection() const
Definition: TGeoMatrix.h:69
virtual void RegisterYourself()
Register the matrix in the current manager, which will become the owner.
Definition: TGeoMatrix.cxx:526
void Print(Option_t *option="") const
print the matrix in 4x4 format
Definition: TGeoMatrix.cxx:486
Bool_t IsIdentity() const
Definition: TGeoMatrix.h:66
Bool_t IsRegistered() const
Definition: TGeoMatrix.h:77
char * GetPointerName() const
Provide a pointer name containing uid.
Definition: TGeoMatrix.cxx:294
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
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Definition: TGeoMedium.cxx:137
void SetMaterial(TGeoMaterial *mat)
Definition: TGeoMedium.h:55
char * GetPointerName() const
Provide a pointer name containing uid.
Definition: TGeoMedium.cxx:127
A node containing local transformation.
Definition: TGeoNode.h:153
void SetMatrix(const TGeoMatrix *matrix)
Matrix setter.
Definition: TGeoNode.cxx:766
virtual TGeoMatrix * GetMatrix() const
Definition: TGeoNode.h:170
Node containing an offset.
Definition: TGeoNode.h:184
void SetFinder(TGeoPatternFinder *finder)
Definition: TGeoNode.h:206
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.
Definition: TGeoShape.cxx:326
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.
Definition: TGeoShape.cxx:699
void CheckShape(Int_t testNo, Int_t nsamples=10000, Option_t *option="")
Test for shape navigation methods.
Definition: TGeoShape.cxx:209
virtual Bool_t IsValidBox() const =0
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:301
static TGeoVolumeAssembly * MakeAssemblyFromVolume(TGeoVolume *vol)
Make a clone of volume VOL but which is an assembly.
virtual Int_t GetCurrentNodeIndex() const
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:317
virtual ~TGeoVolumeAssembly()
Destructor. The assembly is owner of its "shape".
std::mutex fMutex
Thread vector size.
Definition: TGeoVolume.h:319
virtual void AddNode(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat=0, Option_t *option="")
Add a component to the assembly.
Int_t fThreadSize
Thread specific data vector.
Definition: TGeoVolume.h:318
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:252
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:259
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:272
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
Double_t fStep
Definition: TGeoVolume.h:260
virtual void SetLineColor(Color_t lcolor)
Set the line color for all components.
TGeoVolumeMulti * fDivision
Definition: TGeoVolume.h:255
TGeoVolumeMulti()
dummy constructor
virtual void AddNode(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat, Option_t *option="")
Add a new node to the list of nodes.
TGeoShape * GetLastShape() const
Returns the last shape.
Int_t GetNvolumes() const
Definition: TGeoVolume.h:277
virtual void SetMedium(TGeoMedium *medium)
Set medium for a multiple volume.
TObjArray * fVolumes
Definition: TGeoVolume.h:254
virtual void SetLineWidth(Width_t lwidth)
Set the line width for all components.
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition: TGeoVolume.h:47
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.
Definition: TGeoVolume.cxx:970
virtual void cd(Int_t inode) const
Actualize matrix of node indexed <inode>
Definition: TGeoVolume.cxx:922
virtual void ClearThreadData() const
Definition: TGeoVolume.cxx:423
virtual ~TGeoVolume()
Destructor.
Definition: TGeoVolume.cxx:501
virtual void Print(Option_t *option="") const
Print volume info.
Bool_t IsVisContainers() const
Definition: TGeoVolume.h:153
void SetVoxelFinder(TGeoVoxelFinder *finder)
Definition: TGeoVolume.h:228
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.
void RandomRays(Int_t nrays=10000, Double_t startx=0, Double_t starty=0, Double_t startz=0, const char *target_vol=0, Bool_t check_norm=kFALSE)
Random raytracing method.
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.
Definition: TGeoVolume.cxx:752
Int_t fNtotal
Definition: TGeoVolume.h:60
void MakeCopyNodes(const TGeoVolume *other)
make a new list of nodes and copy all nodes of other volume inside
TGeoNode * ReplaceNode(TGeoNode *nodeorig, TGeoShape *newshape=0, TGeoMatrix *newpos=0, TGeoMedium *newmed=0)
Replace an existing daughter with a new volume having the same name but possibly a new shape,...
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:230
void ClearNodes()
Definition: TGeoVolume.h:98
void Raytrace(Bool_t flag=kTRUE)
Draw this volume with current settings and perform raytracing in the pad.
TGeoVolume()
dummy constructor
Definition: TGeoVolume.cxx:447
TGeoMedium * GetMedium() const
Definition: TGeoVolume.h:171
virtual Bool_t IsFolder() const
Return TRUE if volume contains nodes.
Definition: TGeoVolume.cxx:791
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.
Definition: TGeoVolume.cxx:552
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:59
virtual void CreateThreadData(Int_t nthreads)
Definition: TGeoVolume.cxx:431
virtual Int_t GetByteCount() const
get the total size in bytes for this volume
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.
Definition: TGeoVolume.cxx:518
TGeoManager * fGeoManager
Definition: TGeoVolume.h:55
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:54
TGeoMaterial * GetMaterial() const
Definition: TGeoVolume.h:170
virtual Bool_t IsVolumeMulti() const
Definition: TGeoVolume.h:113
TGeoExtension * GrabUserExtension() const
Get a copy of the user extension pointer.
@ kVolumeClone
Definition: TGeoVolume.h:83
@ kVolumeSelected
Definition: TGeoVolume.h:76
@ kVolumeMulti
Definition: TGeoVolume.h:80
@ kVolumeImportNodes
Definition: TGeoVolume.h:79
Int_t CountNodes(Int_t nlevels=1000, Int_t option=0)
Count total number of subnodes starting from this volume, nlevels down.
Definition: TGeoVolume.cxx:701
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:229
Int_t GetNdaughters() const
Definition: TGeoVolume.h:347
Bool_t IsValid() const
Definition: TGeoVolume.h:150
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
void Grab()
Definition: TGeoVolume.h:137
char * GetPointerName() const
Provide a pointer name containing uid.
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...
Definition: TGeoVolume.cxx:567
void SelectVolume(Bool_t clear=kFALSE)
Select this volume as matching an arbitrary criteria.
TObjArray * GetNodes()
Definition: TGeoVolume.h:165
virtual char * GetObjectInfo(Int_t px, Int_t py) const
Get volume info for the browser.
Option_t * GetOption() const
Definition: TGeoVolume.h:183
void ClearShape()
Clear the shape of this volume from the list held by the current manager.
Definition: TGeoVolume.cxx:639
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
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.
Definition: TGeoVolume.cxx:647
void SetNtotal(Int_t ntotal)
Definition: TGeoVolume.h:231
virtual void Paint(Option_t *option="")
paint volume
Bool_t GetOptimalVoxels() const
Returns true if cylindrical voxelization is optimal.
void InvisibleAll(Bool_t flag=kTRUE)
Make volume and each of it daughters (in)visible.
Definition: TGeoVolume.cxx:763
Bool_t IsVisibleDaughters() const
Definition: TGeoVolume.h:152
TString fOption
just a hook for now
Definition: TGeoVolume.h:58
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:212
TGeoPatternFinder * GetFinder() const
Definition: TGeoVolume.h:173
void PrintVoxels() const
Print the voxels for this volume.
TGeoExtension * fUserExtension
Definition: TGeoVolume.h:62
virtual void SetMedium(TGeoMedium *medium)
Definition: TGeoVolume.h:227
TGeoVoxelFinder * GetVoxels() const
Getter for optimization structure.
void SetAttVisibility(Bool_t vis)
Definition: TGeoVolume.h:218
void SetShape(const TGeoShape *shape)
set the shape associated with this volume
static TGeoMedium * DummyMedium()
Definition: TGeoVolume.cxx:439
TObject * fField
pointer to TGeoManager owning this volume
Definition: TGeoVolume.h:57
Int_t GetNumber() const
Definition: TGeoVolume.h:180
void CleanAll()
Clean data of the volume.
Definition: TGeoVolume.cxx:630
Bool_t IsTopVolume() const
True if this is the top volume of the geometry.
Definition: TGeoVolume.cxx:811
TGeoMedium * fMedium
Definition: TGeoVolume.h:51
TGeoShape * GetShape() const
Definition: TGeoVolume.h:186
void InspectMaterial() const
Inspect the material for this volume.
Definition: TGeoVolume.cxx:828
void PrintNodes() const
print nodes
static TGeoMedium * fgDummyMedium
Definition: TGeoVolume.h:52
void RegisterYourself(Option_t *option="")
Register the volume and all materials/media/matrices/shapes to the manager.
virtual void AddNodeOverlap(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat=0, Option_t *option="")
Add a TGeoNode to the list of nodes.
Definition: TGeoVolume.cxx:995
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.
Definition: TGeoVolume.cxx:820
TGeoShape * fShape
Definition: TGeoVolume.h:50
void SetField(TObject *field)
Definition: TGeoVolume.h:216
static TGeoVolume * Import(const char *filename, const char *name="", Option_t *option="")
Import a volume from a file.
Definition: TGeoVolume.cxx:836
Bool_t IsStyleDefault() const
check if the visibility and attributes are the default ones
Definition: TGeoVolume.cxx:799
static void CreateDummyMedium()
Create a dummy medium.
Definition: TGeoVolume.cxx:411
TGeoExtension * fFWExtension
Transient user-defined extension to volumes.
Definition: TGeoVolume.h:63
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:154
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:172
TGeoPatternFinder * fFinder
dummy medium
Definition: TGeoVolume.h:53
Int_t Export(const char *filename, const char *name="", Option_t *option="")
Export this volume to a file.
Definition: TGeoVolume.cxx:886
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.
virtual void AddNode(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat=0, Option_t *option="")
Add a TGeoNode to the list of nodes.
Definition: TGeoVolume.cxx:931
TObjArray * fNodes
Definition: TGeoVolume.h:49
virtual Bool_t IsVisible() const
Definition: TGeoVolume.h:151
Int_t GetNtotal() const
Definition: TGeoVolume.h:167
void InspectShape() const
Definition: TGeoVolume.h:191
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.
Definition: TGeoVolume.cxx:588
void SetOverlappingCandidate(Bool_t flag)
Definition: TGeoVolume.h:213
Bool_t IsOverlappingCandidate() const
Definition: TGeoVolume.h:144
Int_t fRefCount
Definition: TGeoVolume.h:61
void CheckShape(Int_t testNo, Int_t nsamples=10000, Option_t *option="")
Tests for checking the shape navigation algorithms. See TGeoShape::CheckShape()
Definition: TGeoVolume.cxx:622
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.
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2998
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:251
Book space in a file, create I/O buffers, to fill them, (un)compress them.
Definition: TKey.h:24
virtual const char * GetClassName() const
Definition: TKey.h:72
virtual TObject * ReadObj()
To read a TObject* from the file.
Definition: TKey.cxx:729
virtual void Add(TObject *obj)
Definition: TList.h:87
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition: TList.cxx:575
TMap implements an associative array of (key,value) pairs using a THashTable for efficient retrieval ...
Definition: TMap.h:40
void Add(TObject *obj)
This function may not be used (but we need to provide it since it is a pure virtual in TCollection).
Definition: TMap.cxx:53
TObject * GetValue(const char *keyname) const
Returns a pointer to the value associated with keyname as name of the key.
Definition: TMap.cxx:235
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
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
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
An array of TObjects.
Definition: TObjArray.h:37
Int_t GetEntriesFast() const
Definition: TObjArray.h:64
virtual void AddAtAndExpand(TObject *obj, Int_t idx)
Add object at position idx.
Definition: TObjArray.cxx:234
void Add(TObject *obj)
Definition: TObjArray.h:74
virtual void Compress()
Remove empty slots from array.
Definition: TObjArray.cxx:333
Int_t GetEntries() const
Return the number of objects in array (i.e.
Definition: TObjArray.cxx:522
virtual void Clear(Option_t *option="")
Remove all objects from the array.
Definition: TObjArray.cxx:320
virtual void Delete(Option_t *option="")
Remove all objects from the array AND delete all heap based objects.
Definition: TObjArray.cxx:355
virtual TObject * FindObject(const char *name) const
Find an object in this collection using its name.
Definition: TObjArray.cxx:414
virtual TObject * Remove(TObject *obj)
Remove object from array.
Definition: TObjArray.cxx:718
virtual void AddAt(TObject *obj, Int_t idx)
Add object at position ids.
Definition: TObjArray.cxx:253
virtual TObject * RemoveAt(Int_t idx)
Remove object at index idx.
Definition: TObjArray.cxx:693
TObject * At(Int_t idx) const
Definition: TObjArray.h:166
Mother of all ROOT objects.
Definition: TObject.h:37
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition: TObject.cxx:785
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition: TObject.h:172
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:128
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:866
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Definition: TObject.cxx:664
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:694
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:443
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition: TObject.cxx:908
virtual void SetUniqueID(UInt_t uid)
Set the unique object id.
Definition: TObject.cxx:705
void ResetBit(UInt_t f)
Definition: TObject.h:171
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition: TObject.cxx:854
Sequenceable collection abstract base class.
Basic string class.
Definition: TString.h:131
Ssiz_t Length() const
Definition: TString.h:405
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1125
TSubString Strip(EStripType s=kTrailing, char c=' ') const
Return a substring of self stripped at beginning and/or end.
Definition: TString.cxx:1106
const char * Data() const
Definition: TString.h:364
Ssiz_t Capacity() const
Definition: TString.h:352
Bool_t IsNull() const
Definition: TString.h:402
TString & Remove(Ssiz_t pos)
Definition: TString.h:668
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:2311
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:619
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:634
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 TH2F * LegoPlot(Int_t ntheta=60, Double_t themin=0., Double_t themax=180., Int_t nphi=90, Double_t phimin=0., Double_t phimax=360., Double_t rmin=0., Double_t rmax=9999999, Option_t *option="")=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
double dist(Rotation3D const &r1, Rotation3D const &r2)
Definition: 3DDistances.cxx:48
constexpr Double_t E()
Base of natural log:
Definition: TMath.h:97
ThreadData_t()
index of next node to be entered