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