Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
TGeoOverlap.cxx
Go to the documentation of this file.
1 // @(#)root/geom:$Id$
2 // Author: Andrei Gheata 09-02-03
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 #include "TVirtualPad.h"
13 #include "TMath.h"
14 #include "TNamed.h"
15 #include "TBrowser.h"
16 #include "TGeoManager.h"
17 #include "TGeoVolume.h"
18 #include "TGeoNode.h"
19 #include "TGeoBBox.h"
20 #include "TRandom3.h"
21 #include "TPolyMarker3D.h"
22 #include "TVirtualGeoPainter.h"
23 
24 #include "TGeoOverlap.h"
25 
26 ClassImp(TGeoOverlap);
27 
28 /** \class TGeoOverlap
29 \ingroup Geometry_classes
30 
31 Base class describing geometry overlaps. Overlaps apply
32 to the nodes contained inside a volume. These should not overlap to
33 each other nor extrude the shape of their mother volume.
34 */
35 
36 ////////////////////////////////////////////////////////////////////////////////
37 /// Default ctor.
38 
39 TGeoOverlap::TGeoOverlap()
40 {
41  fOverlap = 0;
42  fVolume1 = 0;
43  fVolume2 = 0;
44  fMatrix1 = 0;
45  fMatrix2 = 0;
46  fMarker = 0;
47 }
48 
49 ////////////////////////////////////////////////////////////////////////////////
50 /// Creates a named overlap belonging to volume VOL and having the size OVLP.
51 
52 TGeoOverlap::TGeoOverlap(const char *name, TGeoVolume *vol1, TGeoVolume *vol2,
53  const TGeoMatrix *matrix1, const TGeoMatrix *matrix2,
54  Bool_t isovlp, Double_t ovlp)
55  :TNamed("",name)
56 {
57  fOverlap = ovlp;
58  fVolume1 = vol1;
59  fVolume2 = vol2;
60  fMatrix1 = new TGeoHMatrix();
61  *fMatrix1 = matrix1;
62  fMatrix2 = new TGeoHMatrix();
63  *fMatrix2 = matrix2;
64  fMarker = new TPolyMarker3D();
65  fMarker->SetMarkerColor(2);
66  SetIsOverlap(isovlp);
67  fMarker->SetMarkerStyle(6);
68 // fMarker->SetMarkerSize(0.5);
69 }
70 
71 ////////////////////////////////////////////////////////////////////////////////
72 /// Destructor.
73 
74 TGeoOverlap::~TGeoOverlap()
75 {
76  if (fMarker) delete fMarker;
77  if (fMatrix1) delete fMatrix1;
78  if (fMatrix2) delete fMatrix2;
79 }
80 
81 ////////////////////////////////////////////////////////////////////////////////
82 /// Define double-click action
83 
84 void TGeoOverlap::Browse(TBrowser *b)
85 {
86  if (!b) return;
87  Draw();
88 }
89 
90 ////////////////////////////////////////////////////////////////////////////////
91 /// Method to compare this overlap with another. Returns :
92 /// -1 - this is smaller than OBJ
93 /// 0 - equal
94 /// 1 - greater
95 
96 Int_t TGeoOverlap::Compare(const TObject *obj) const
97 {
98  TGeoOverlap *other = 0;
99  other = (TGeoOverlap*)obj;
100  if (!other) {
101  Error("Compare", "other object is not TGeoOverlap");
102  return 0;
103  }
104  if (IsExtrusion()) {
105  if (other->IsExtrusion()) return (fOverlap<=other->GetOverlap())?1:-1;
106  return -1;
107  } else {
108  if (other->IsExtrusion()) return 1;
109  return (fOverlap<=other->GetOverlap())?1:-1;
110  }
111 }
112 
113 ////////////////////////////////////////////////////////////////////////////////
114 /// Distance to primitive for an overlap.
115 
116 Int_t TGeoOverlap::DistancetoPrimitive(Int_t px, Int_t py)
117 {
118  return fVolume1->GetGeoManager()->GetGeomPainter()->DistanceToPrimitiveVol(fVolume1, px, py);
119 }
120 
121 ////////////////////////////////////////////////////////////////////////////////
122 /// Draw the overlap. One daughter will be blue, the other green,
123 /// extruding points red.
124 
125 void TGeoOverlap::Draw(Option_t *option)
126 {
127  fVolume1->GetGeoManager()->GetGeomPainter()->DrawOverlap(this, option);
128  PrintInfo();
129 }
130 
131 ////////////////////////////////////////////////////////////////////////////////
132 /// Event interception.
133 
134 void TGeoOverlap::ExecuteEvent(Int_t event, Int_t px, Int_t py)
135 {
136  fVolume1->GetGeoManager()->GetGeomPainter()->ExecuteVolumeEvent(fVolume1, event, px, py);
137 }
138 
139 ////////////////////////////////////////////////////////////////////////////////
140 /// Paint the overlap.
141 
142 void TGeoOverlap::Paint(Option_t *option)
143 {
144  fVolume1->GetGeoManager()->GetGeomPainter()->PaintOverlap(this, option);
145 }
146 
147 ////////////////////////////////////////////////////////////////////////////////
148 /// Print detailed info.
149 
150 void TGeoOverlap::Print(Option_t *) const
151 {
152  PrintInfo();
153  printf(" - first volume: %s at position:\n", fVolume1->GetName());
154  fMatrix1->Print();
155  fVolume1->InspectShape();
156  printf(" - second volume: %s at position:\n", fVolume2->GetName());
157  fMatrix2->Print();
158  fVolume2->InspectShape();
159 }
160 
161 ////////////////////////////////////////////////////////////////////////////////
162 /// Print some info.
163 
164 void TGeoOverlap::PrintInfo() const
165 {
166  printf(" = Overlap %s: %s ovlp=%g\n", GetName(), GetTitle(),fOverlap);
167 }
168 
169 ////////////////////////////////////////////////////////////////////////////////
170 /// Set next overlapping point.
171 
172 void TGeoOverlap::SetNextPoint(Double_t x, Double_t y, Double_t z)
173 {
174  fMarker->SetNextPoint(x,y,z);
175 }
176 
177 ////////////////////////////////////////////////////////////////////////////////
178 /// Draw overlap and sample with random points the overlapping region.
179 
180 void TGeoOverlap::SampleOverlap(Int_t npoints)
181 {
182  Draw();
183  // Select bounding box of the second volume (may extrude first)
184  TPolyMarker3D *marker = 0;
185  TGeoBBox *box = (TGeoBBox*)fVolume2->GetShape();
186  Double_t dx = box->GetDX();
187  Double_t dy = box->GetDY();
188  Double_t dz = box->GetDZ();
189  Double_t pt[3];
190  Double_t master[3];
191  const Double_t *orig = box->GetOrigin();
192  Int_t ipoint = 0;
193  Int_t itry = 0;
194  Int_t iovlp = 0;
195  while (ipoint < npoints) {
196  // Shoot randomly in the bounding box.
197  pt[0] = orig[0] - dx + 2.*dx*gRandom->Rndm();
198  pt[1] = orig[1] - dy + 2.*dy*gRandom->Rndm();
199  pt[2] = orig[2] - dz + 2.*dz*gRandom->Rndm();
200  if (!fVolume2->Contains(pt)) {
201  itry++;
202  if (itry>10000 && !ipoint) {
203  Error("SampleOverlap", "No point inside volume!!! - aborting");
204  break;
205  }
206  continue;
207  }
208  ipoint++;
209  // Check if the point is inside the first volume
210  fMatrix2->LocalToMaster(pt, master);
211  fMatrix1->MasterToLocal(master, pt);
212  Bool_t in = fVolume1->Contains(pt);
213  if (IsOverlap() && !in) continue;
214  if (!IsOverlap() && in) continue;
215  // The point is in the overlapping region.
216  iovlp++;
217  if (!marker) {
218  marker = new TPolyMarker3D();
219  marker->SetMarkerColor(kRed);
220  }
221  marker->SetNextPoint(master[0], master[1], master[2]);
222  }
223  if (!iovlp) return;
224  marker->Draw("SAME");
225  gPad->Modified();
226  gPad->Update();
227  Double_t capacity = fVolume1->GetShape()->Capacity();
228  capacity *= Double_t(iovlp)/Double_t(npoints);
229  Double_t err = 1./TMath::Sqrt(Double_t(iovlp));
230  Info("SampleOverlap", "#Overlap %s has %g +/- %g [cm3]",
231  GetName(), capacity, err*capacity);
232 }
233 
234 ////////////////////////////////////////////////////////////////////////////////
235 /// Get 3D size of this.
236 
237 void TGeoOverlap::Sizeof3D() const
238 {
239  fVolume1->GetShape()->Sizeof3D();
240  fVolume2->GetShape()->Sizeof3D();
241 }
242 
243 ////////////////////////////////////////////////////////////////////////////////
244 /// Validate this overlap.
245 
246 void TGeoOverlap::Validate() const
247 {
248  Double_t point[3];
249  Double_t local[3];
250  Double_t safe1,safe2;
251  Int_t npoints = fMarker->GetN();
252  for (Int_t i=0; i<npoints; i++) {
253  fMarker->GetPoint(i, point[0], point[1], point[2]);
254  if (IsExtrusion()) {
255  fMatrix1->MasterToLocal(point,local);
256  safe1 = fVolume1->GetShape()->Safety(local, kFALSE);
257  printf("point %d: safe1=%f\n", i, safe1);
258  } else {
259  fMatrix1->MasterToLocal(point,local);
260  safe1 = fVolume1->GetShape()->Safety(local, kTRUE);
261  fMatrix2->MasterToLocal(point,local);
262  safe2 = fVolume2->GetShape()->Safety(local, kTRUE);
263  printf("point %d: safe1=%f safe2=%f\n", i, safe1,safe2);
264  }
265  }
266 }
267 
268