Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoPara.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Andrei Gheata 31/01/02
3// TGeoPara::Contains() 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 TGeoPara
14\ingroup Shapes_classes
15\brief Parallelepiped class.
16
17A parallelepiped is a shape having 3 pairs of parallel faces out of
18which one is parallel with the XY plane (Z faces). All faces are
19parallelograms in the general case. The Z faces have 2 edges parallel
20with the X-axis.
21
22The shape has the center in the origin and it is defined by:
23
24 - `dX, dY, dZ:` half-lengths of the projections of the edges on X, Y
25 and Z. The lower Z face is positioned at `-dZ`, while the upper at `+dZ`.
26 - `alpha:` angle between the segment defined by the centers of the
27 X-parallel edges and Y axis `[-90,90]` in degrees
28 - `theta:` theta angle of the segment defined by the centers of the Z faces;
29 - `phi:` phi angle of the same segment
30
31~~~ {.cpp}
32TGeoPara(dX,dY,dZ,alpha,theta,phi);
33~~~
34
35A box is a particular parallelepiped having the parameters:
36`(dX,dY,dZ,0.,0.,0.)`.
37
38Begin_Macro
39{
40 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
41 new TGeoManager("para", "poza1");
42 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
43 TGeoMedium *med = new TGeoMedium("MED",1,mat);
44 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
45 gGeoManager->SetTopVolume(top);
46 TGeoVolume *vol = gGeoManager->MakePara("PARA",med, 20,30,40,30,15,30);
47 vol->SetLineWidth(2);
48 top->AddNode(vol,1);
49 gGeoManager->CloseGeometry();
50 gGeoManager->SetNsegments(80);
51 top->Draw();
52 TView *view = gPad->GetView();
53 view->ShowAxis();
54}
55End_Macro
56
57*/
58
59#include <iostream>
60
61#include "TGeoManager.h"
62#include "TGeoMatrix.h"
63#include "TGeoVolume.h"
64#include "TGeoPara.h"
65#include "TMath.h"
66
68
69////////////////////////////////////////////////////////////////////////////////
70/// Default constructor
71
73{
75 fX = fY = fZ = 0;
76 fAlpha = 0;
77 fTheta = 0;
78 fPhi = 0;
79 fTxy = 0;
80 fTxz = 0;
81 fTyz = 0;
82}
83
84////////////////////////////////////////////////////////////////////////////////
85/// Default constructor specifying minimum and maximum radius
86
88 Double_t theta, Double_t phi)
89 :TGeoBBox(0, 0, 0)
90{
92 fX = dx;
93 fY = dy;
94 fZ = dz;
95 fAlpha = alpha;
96 fTheta = theta;
97 fPhi = phi;
99 Double_t tth = TMath::Tan(theta*TMath::DegToRad());
100 Double_t ph = phi*TMath::DegToRad();
101 fTxz = tth*TMath::Cos(ph);
102 fTyz = tth*TMath::Sin(ph);
103 if ((fX<0) || (fY<0) || (fZ<0)) {
104// printf("para : %f %f %f\n", fX, fY, fZ);
106 }
107 else ComputeBBox();
108}
109
110////////////////////////////////////////////////////////////////////////////////
111/// Default constructor specifying minimum and maximum radius
112
114 Double_t theta, Double_t phi)
115 :TGeoBBox(name, 0, 0, 0)
116{
118 fX = dx;
119 fY = dy;
120 fZ = dz;
121 fAlpha = alpha;
122 fTheta = theta;
123 fPhi = phi;
125 Double_t tth = TMath::Tan(theta*TMath::DegToRad());
126 Double_t ph = phi*TMath::DegToRad();
127 fTxz = tth*TMath::Cos(ph);
128 fTyz = tth*TMath::Sin(ph);
129 if ((fX<0) || (fY<0) || (fZ<0)) {
130// printf("para : %f %f %f\n", fX, fY, fZ);
132 }
133 else ComputeBBox();
134}
135
136////////////////////////////////////////////////////////////////////////////////
137/// Default constructor
138/// - param[0] = dx
139/// - param[1] = dy
140/// - param[2] = dz
141/// - param[3] = alpha
142/// - param[4] = theta
143/// - param[5] = phi
144
146 :TGeoBBox(0, 0, 0)
147{
149 SetDimensions(param);
150 if ((fX<0) || (fY<0) || (fZ<0)) SetShapeBit(kGeoRunTimeShape);
151 else ComputeBBox();
152}
153
154////////////////////////////////////////////////////////////////////////////////
155/// destructor
156
158{
159}
160
161////////////////////////////////////////////////////////////////////////////////
162/// Computes capacity of the shape in [length^3]
163
165{
166 Double_t capacity = 8.*fX*fY*fZ;
167 return capacity;
168}
169
170////////////////////////////////////////////////////////////////////////////////
171/// compute bounding box
172
174{
177 Double_t dz = fZ;
178 TGeoBBox::SetBoxDimensions(dx, dy, dz);
179 memset(fOrigin, 0, 3*sizeof(Double_t));
180}
181
182////////////////////////////////////////////////////////////////////////////////
183/// Compute normal to closest surface from POINT.
184
185void TGeoPara::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
186{
187 Double_t saf[3];
188 // distance from point to higher Z face
189 saf[0] = TMath::Abs(fZ-TMath::Abs(point[2])); // Z
190
191 Double_t yt = point[1]-fTyz*point[2];
192 saf[1] = TMath::Abs(fY-TMath::Abs(yt)); // Y
193 // cos of angle YZ
194 Double_t cty = 1.0/TMath::Sqrt(1.0+fTyz*fTyz);
195
196 Double_t xt = point[0]-fTxz*point[2]-fTxy*yt;
197 saf[2] = TMath::Abs(fX-TMath::Abs(xt)); // X
198 // cos of angle XZ
199 Double_t ctx = 1.0/TMath::Sqrt(1.0+fTxy*fTxy+fTxz*fTxz);
200 saf[2] *= ctx;
201 saf[1] *= cty;
202 Int_t i = TMath::LocMin(3,saf);
203 switch (i) {
204 case 0:
205 norm[0] = norm[1] = 0;
206 norm[2] = TMath::Sign(1.,dir[2]);
207 return;
208 case 1:
209 norm[0] = 0;
210 norm[1] = cty;
211 norm[2] = - fTyz*cty;
212 break;
213 case 2:
216 norm[2] = -TMath::Sin(fTheta*TMath::DegToRad());
217 }
218 if (norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2]<0) {
219 norm[0] = -norm[0];
220 norm[1] = -norm[1];
221 norm[2] = -norm[2];
222 }
223}
224
225////////////////////////////////////////////////////////////////////////////////
226/// test if point is inside this sphere
227/// test Z range
228
230{
231 if (TMath::Abs(point[2]) > fZ) return kFALSE;
232 // check X and Y
233 Double_t yt=point[1]-fTyz*point[2];
234 if (TMath::Abs(yt) > fY) return kFALSE;
235 Double_t xt=point[0]-fTxz*point[2]-fTxy*yt;
236 if (TMath::Abs(xt) > fX) return kFALSE;
237 return kTRUE;
238}
239
240////////////////////////////////////////////////////////////////////////////////
241/// compute distance from inside point to surface of the para
242/// Boundary safe algorithm.
243
244Double_t TGeoPara::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
245{
246 if (iact<3 && safe) {
247 // compute safety
248 *safe = Safety(point, kTRUE);
249 if (iact==0) return TGeoShape::Big();
250 if (iact==1 && step<*safe) return TGeoShape::Big();
251 }
252 Double_t saf[2];
253 Double_t snxt = TGeoShape::Big();
254 Double_t s;
255 saf[0] = fZ+point[2];
256 saf[1] = fZ-point[2];
257 if (!TGeoShape::IsSameWithinTolerance(dir[2],0)) {
258 s = (dir[2]>0)?(saf[1]/dir[2]):(-saf[0]/dir[2]);
259 if (s<0) return 0.0;
260 if (s<snxt) snxt = s;
261 }
262 // distance from point to center axis on Y
263 Double_t yt = point[1]-fTyz*point[2];
264 saf[0] = fY+yt;
265 saf[1] = fY-yt;
266 Double_t dy = dir[1]-fTyz*dir[2];
268 s = (dy>0)?(saf[1]/dy):(-saf[0]/dy);
269 if (s<0) return 0.0;
270 if (s<snxt) snxt = s;
271 }
272 // distance from point to center axis on X
273 Double_t xt = point[0]-fTxz*point[2]-fTxy*yt;
274 saf[0] = fX+xt;
275 saf[1] = fX-xt;
276 Double_t dx = dir[0]-fTxz*dir[2]-fTxy*dy;
278 s = (dx>0)?(saf[1]/dx):(-saf[0]/dx);
279 if (s<0) return 0.0;
280 if (s<snxt) snxt = s;
281 }
282 return snxt;
283}
284
285////////////////////////////////////////////////////////////////////////////////
286/// compute distance from inside point to surface of the para
287
288Double_t TGeoPara::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
289{
290 if (iact<3 && safe) {
291 // compute safe distance
292 *safe = Safety(point, kFALSE);
293 if (iact==0) return TGeoShape::Big();
294 if (iact==1 && step<*safe) return TGeoShape::Big();
295 }
296 Bool_t in = kTRUE;
297 Double_t safz;
298 safz = TMath::Abs(point[2])-fZ;
299 if (safz>0) {
300 // outside Z
301 if (point[2]*dir[2]>=0) return TGeoShape::Big();
302 in = kFALSE;
303 }
304 Double_t yt=point[1]-fTyz*point[2];
305 Double_t safy = TMath::Abs(yt)-fY;
306 Double_t dy=dir[1]-fTyz*dir[2];
307 if (safy>0) {
308 if (yt*dy>=0) return TGeoShape::Big();
309 in = kFALSE;
310 }
311 Double_t xt=point[0]-fTxy*yt-fTxz*point[2];
312 Double_t safx = TMath::Abs(xt)-fX;
313 Double_t dx=dir[0]-fTxy*dy-fTxz*dir[2];
314 if (safx>0) {
315 if (xt*dx>=0) return TGeoShape::Big();
316 in = kFALSE;
317 }
318 // protection in case point is actually inside
319 if (in) {
320 if (safz>safx && safz>safy) {
321 if (point[2]*dir[2]>0) return TGeoShape::Big();
322 return 0.0;
323 }
324 if (safx>safy) {
325 if (xt*dx>0) return TGeoShape::Big();
326 return 0.0;
327 }
328 if (yt*dy>0) return TGeoShape::Big();
329 return 0.0;
330 }
331 Double_t xnew,ynew,znew;
332 if (safz>0) {
333 Double_t snxt = safz/TMath::Abs(dir[2]);
334 xnew = point[0]+snxt*dir[0];
335 ynew = point[1]+snxt*dir[1];
336 znew = (point[2]>0)?fZ:(-fZ);
337 Double_t ytn = ynew-fTyz*znew;
338 if (TMath::Abs(ytn)<=fY) {
339 Double_t xtn = xnew-fTxy*ytn-fTxz*znew;
340 if (TMath::Abs(xtn)<=fX) return snxt;
341 }
342 }
343 if (safy>0) {
344 Double_t snxt = safy/TMath::Abs(dy);
345 znew = point[2]+snxt*dir[2];
346 if (TMath::Abs(znew)<=fZ) {
347 Double_t ytn = (yt>0)?fY:(-fY);
348 xnew = point[0]+snxt*dir[0];
349 Double_t xtn = xnew-fTxy*ytn-fTxz*znew;
350 if (TMath::Abs(xtn)<=fX) return snxt;
351 }
352 }
353 if (safx>0) {
354 Double_t snxt = safx/TMath::Abs(dx);
355 znew = point[2]+snxt*dir[2];
356 if (TMath::Abs(znew)<=fZ) {
357 ynew = point[1]+snxt*dir[1];
358 Double_t ytn = ynew-fTyz*znew;
359 if (TMath::Abs(ytn)<=fY) return snxt;
360 }
361 }
362 return TGeoShape::Big();
363}
364
365////////////////////////////////////////////////////////////////////////////////
366/// Divide this parallelepiped shape belonging to volume "voldiv" into ndiv equal volumes
367/// called divname, from start position with the given step. Returns pointer
368/// to created division cell volume. In case a wrong division axis is supplied,
369/// returns pointer to volume to be divided.
370
371TGeoVolume *TGeoPara::Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv,
372 Double_t start, Double_t step)
373{
374 TGeoShape *shape; //--- shape to be created
375 TGeoVolume *vol; //--- division volume to be created
376 TGeoVolumeMulti *vmulti; //--- generic divided volume
377 TGeoPatternFinder *finder; //--- finder to be attached
378 TString opt = ""; //--- option to be attached
379 Double_t end=start+ndiv*step;
380 switch (iaxis) {
381 case 1: //--- divide on X
382 shape = new TGeoPara(step/2, fY, fZ,fAlpha,fTheta, fPhi);
383 finder = new TGeoPatternParaX(voldiv, ndiv, start, end);
384 opt = "X";
385 break;
386 case 2: //--- divide on Y
387 shape = new TGeoPara(fX, step/2, fZ, fAlpha, fTheta, fPhi);
388 finder = new TGeoPatternParaY(voldiv, ndiv, start, end);
389 opt = "Y";
390 break;
391 case 3: //--- divide on Z
392 shape = new TGeoPara(fX, fY, step/2, fAlpha, fTheta, fPhi);
393 finder = new TGeoPatternParaZ(voldiv, ndiv, start, end);
394 opt = "Z";
395 break;
396 default:
397 Error("Divide", "Wrong axis type for division");
398 return 0;
399 }
400 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
401 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
402 vmulti->AddVolume(vol);
403 voldiv->SetFinder(finder);
404 finder->SetDivIndex(voldiv->GetNdaughters());
405 for (Int_t ic=0; ic<ndiv; ic++) {
406 voldiv->AddNodeOffset(vol, ic, start+step/2.+ic*step, opt.Data());
407 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
408 }
409 return vmulti;
410}
411
412////////////////////////////////////////////////////////////////////////////////
413/// Get range of shape for a given axis.
414
416{
417 xlo = 0;
418 xhi = 0;
419 Double_t dx = 0;
420 switch (iaxis) {
421 case 1:
422 xlo = -fX;
423 xhi = fX;
424 dx = xhi-xlo;
425 return dx;
426 case 2:
427 xlo = -fY;
428 xhi = fY;
429 dx = xhi-xlo;
430 return dx;
431 case 3:
432 xlo = -fZ;
433 xhi = fZ;
434 dx = xhi-xlo;
435 return dx;
436 }
437 return dx;
438}
439
440////////////////////////////////////////////////////////////////////////////////
441/// Fill vector param[4] with the bounding cylinder parameters. The order
442/// is the following : Rmin, Rmax, Phi1, Phi2
443
445{
447}
448
449////////////////////////////////////////////////////////////////////////////////
450/// Fills real parameters of a positioned box inside this. Returns 0 if successful.
451
452Int_t TGeoPara::GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const
453{
454 dx=dy=dz=0;
455 if (mat->IsRotation()) {
456 Error("GetFittingBox", "cannot handle parametrized rotated volumes");
457 return 1; // ### rotation not accepted ###
458 }
459 //--> translate the origin of the parametrized box to the frame of this box.
460 Double_t origin[3];
461 mat->LocalToMaster(parambox->GetOrigin(), origin);
462 if (!Contains(origin)) {
463 Error("GetFittingBox", "wrong matrix - parametrized box is outside this");
464 return 1; // ### wrong matrix ###
465 }
466 //--> now we have to get the valid range for all parametrized axis
467 Double_t dd[3];
468 dd[0] = parambox->GetDX();
469 dd[1] = parambox->GetDY();
470 dd[2] = parambox->GetDZ();
471 //-> check if Z range is fixed
472 if (dd[2]<0) {
473 dd[2] = TMath::Min(origin[2]+fZ, fZ-origin[2]);
474 if (dd[2]<0) {
475 Error("GetFittingBox", "wrong matrix");
476 return 1;
477 }
478 }
479 if (dd[0]>=0 && dd[1]>=0) {
480 dx = dd[0];
481 dy = dd[1];
482 dz = dd[2];
483 return 0;
484 }
485 //-> check now range at Z = origin[2] +/- dd[2]
486 Double_t upper[8];
487 Double_t lower[8];
488 Double_t z=origin[2]-dd[2];
489 lower[0]=z*fTxz-fTxy*fY-fX;
490 lower[1]=-fY+z*fTyz;
491 lower[2]=z*fTxz+fTxy*fY-fX;
492 lower[3]=fY+z*fTyz;
493 lower[4]=z*fTxz+fTxy*fY+fX;
494 lower[5]=fY+z*fTyz;
495 lower[6]=z*fTxz-fTxy*fY+fX;
496 lower[7]=-fY+z*fTyz;
497 z=origin[2]+dd[2];
498 upper[0]=z*fTxz-fTxy*fY-fX;
499 upper[1]=-fY+z*fTyz;
500 upper[2]=z*fTxz+fTxy*fY-fX;
501 upper[3]=fY+z*fTyz;
502 upper[4]=z*fTxz+fTxy*fY+fX;
503 upper[5]=fY+z*fTyz;
504 upper[6]=z*fTxz-fTxy*fY+fX;
505 upper[7]=-fY+z*fTyz;
506
507 for (Int_t iaxis=0; iaxis<2; iaxis++) {
508 if (dd[iaxis]>=0) continue;
509 Double_t ddmin = TGeoShape::Big();
510 for (Int_t ivert=0; ivert<4; ivert++) {
511 ddmin = TMath::Min(ddmin, TMath::Abs(origin[iaxis]-lower[2*ivert+iaxis]));
512 ddmin = TMath::Min(ddmin, TMath::Abs(origin[iaxis]-upper[2*ivert+iaxis]));
513 }
514 dd[iaxis] = ddmin;
515 }
516 dx = dd[0];
517 dy = dd[1];
518 dz = dd[2];
519 return 0;
520}
521
522////////////////////////////////////////////////////////////////////////////////
523/// in case shape has some negative parameters, these has to be computed
524/// in order to fit the mother
525
527{
528 if (!TestShapeBit(kGeoRunTimeShape)) return 0;
529 if (!mother->TestShapeBit(kGeoPara)) {
530 Error("GetMakeRuntimeShape", "invalid mother");
531 return 0;
532 }
533 Double_t dx, dy, dz;
534 if (fX<0) dx=((TGeoPara*)mother)->GetX();
535 else dx=fX;
536 if (fY<0) dy=((TGeoPara*)mother)->GetY();
537 else dy=fY;
538 if (fZ<0) dz=((TGeoPara*)mother)->GetZ();
539 else dz=fZ;
540 return (new TGeoPara(dx, dy, dz, fAlpha, fTheta, fPhi));
541}
542
543////////////////////////////////////////////////////////////////////////////////
544/// print shape parameters
545
547{
548 printf("*** Shape %s: TGeoPara ***\n", GetName());
549 printf(" dX = %11.5f\n", fX);
550 printf(" dY = %11.5f\n", fY);
551 printf(" dZ = %11.5f\n", fZ);
552 printf(" alpha = %11.5f\n", fAlpha);
553 printf(" theta = %11.5f\n", fTheta);
554 printf(" phi = %11.5f\n", fPhi);
555 printf(" Bounding box:\n");
557}
558
559////////////////////////////////////////////////////////////////////////////////
560/// computes the closest distance from given point to this shape, according
561/// to option. The matching point on the shape is stored in spoint.
562
564{
565 Double_t saf[3];
566 // distance from point to higher Z face
567 saf[0] = fZ-TMath::Abs(point[2]); // Z
568
569 Double_t yt = point[1]-fTyz*point[2];
570 saf[1] = fY-TMath::Abs(yt); // Y
571 // cos of angle YZ
572 Double_t cty = 1.0/TMath::Sqrt(1.0+fTyz*fTyz);
573
574 Double_t xt = point[0]-fTxz*point[2]-fTxy*yt;
575 saf[2] = fX-TMath::Abs(xt); // X
576 // cos of angle XZ
577 Double_t ctx = 1.0/TMath::Sqrt(1.0+fTxy*fTxy+fTxz*fTxz);
578 saf[2] *= ctx;
579 saf[1] *= cty;
580 if (in) return saf[TMath::LocMin(3,saf)];
581 for (Int_t i=0; i<3; i++) saf[i]=-saf[i];
582 return saf[TMath::LocMax(3,saf)];
583}
584
585////////////////////////////////////////////////////////////////////////////////
586/// Save a primitive as a C++ statement(s) on output stream "out".
587
588void TGeoPara::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
589{
591 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
592 out << " dx = " << fX << ";" << std::endl;
593 out << " dy = " << fY << ";" << std::endl;
594 out << " dz = " << fZ << ";" << std::endl;
595 out << " alpha = " << fAlpha<< ";" << std::endl;
596 out << " theta = " << fTheta << ";" << std::endl;
597 out << " phi = " << fPhi << ";" << std::endl;
598 out << " TGeoShape *" << GetPointerName() << " = new TGeoPara(\"" << GetName() << "\",dx,dy,dz,alpha,theta,phi);" << std::endl;
600}
601
602////////////////////////////////////////////////////////////////////////////////
603/// Set dimensions starting from an array.
604
606{
607 fX = param[0];
608 fY = param[1];
609 fZ = param[2];
610 fAlpha = param[3];
611 fTheta = param[4];
612 fPhi = param[5];
613 fTxy = TMath::Tan(param[3]*TMath::DegToRad());
614 Double_t tth = TMath::Tan(param[4]*TMath::DegToRad());
615 Double_t ph = param[5]*TMath::DegToRad();
616 fTxz = tth*TMath::Cos(ph);
617 fTyz = tth*TMath::Sin(ph);
618}
619
620////////////////////////////////////////////////////////////////////////////////
621/// Create PARA mesh points
622
624{
625 if (!points) return;
626 Double_t txy = fTxy;
627 Double_t txz = fTxz;
628 Double_t tyz = fTyz;
629 *points++ = -fZ*txz-txy*fY-fX; *points++ = -fY-fZ*tyz; *points++ = -fZ;
630 *points++ = -fZ*txz+txy*fY-fX; *points++ = +fY-fZ*tyz; *points++ = -fZ;
631 *points++ = -fZ*txz+txy*fY+fX; *points++ = +fY-fZ*tyz; *points++ = -fZ;
632 *points++ = -fZ*txz-txy*fY+fX; *points++ = -fY-fZ*tyz; *points++ = -fZ;
633 *points++ = +fZ*txz-txy*fY-fX; *points++ = -fY+fZ*tyz; *points++ = +fZ;
634 *points++ = +fZ*txz+txy*fY-fX; *points++ = +fY+fZ*tyz; *points++ = +fZ;
635 *points++ = +fZ*txz+txy*fY+fX; *points++ = +fY+fZ*tyz; *points++ = +fZ;
636 *points++ = +fZ*txz-txy*fY+fX; *points++ = -fY+fZ*tyz; *points++ = +fZ;
637}
638
639////////////////////////////////////////////////////////////////////////////////
640/// create sphere mesh points
641
643{
644 if (!points) return;
645 Double_t txy = fTxy;
646 Double_t txz = fTxz;
647 Double_t tyz = fTyz;
648 *points++ = -fZ*txz-txy*fY-fX; *points++ = -fY-fZ*tyz; *points++ = -fZ;
649 *points++ = -fZ*txz+txy*fY-fX; *points++ = +fY-fZ*tyz; *points++ = -fZ;
650 *points++ = -fZ*txz+txy*fY+fX; *points++ = +fY-fZ*tyz; *points++ = -fZ;
651 *points++ = -fZ*txz-txy*fY+fX; *points++ = -fY-fZ*tyz; *points++ = -fZ;
652 *points++ = +fZ*txz-txy*fY-fX; *points++ = -fY+fZ*tyz; *points++ = +fZ;
653 *points++ = +fZ*txz+txy*fY-fX; *points++ = +fY+fZ*tyz; *points++ = +fZ;
654 *points++ = +fZ*txz+txy*fY+fX; *points++ = +fY+fZ*tyz; *points++ = +fZ;
655 *points++ = +fZ*txz-txy*fY+fX; *points++ = -fY+fZ*tyz; *points++ = +fZ;
656}
657
658////////////////////////////////////////////////////////////////////////////////
659/// fill size of this 3-D object
660
662{
664}
665
666////////////////////////////////////////////////////////////////////////////////
667/// Check the inside status for each of the points in the array.
668/// Input: Array of point coordinates + vector size
669/// Output: Array of Booleans for the inside of each point
670
671void TGeoPara::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
672{
673 for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
674}
675
676////////////////////////////////////////////////////////////////////////////////
677/// Compute the normal for an array o points so that norm.dot.dir is positive
678/// Input: Arrays of point coordinates and directions + vector size
679/// Output: Array of normal directions
680
681void TGeoPara::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
682{
683 for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
684}
685
686////////////////////////////////////////////////////////////////////////////////
687/// Compute distance from array of input points having directions specified by dirs. Store output in dists
688
689void TGeoPara::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
690{
691 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
692}
693
694////////////////////////////////////////////////////////////////////////////////
695/// Compute distance from array of input points having directions specified by dirs. Store output in dists
696
697void TGeoPara::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
698{
699 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
700}
701
702////////////////////////////////////////////////////////////////////////////////
703/// Compute safe distance from each of the points in the input array.
704/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
705/// Output: Safety values
706
707void TGeoPara::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
708{
709 for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
710}
const Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
float Float_t
Definition RtypesCore.h:57
const Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:364
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
point * points
Definition X3DBuffer.c:22
Box class.
Definition TGeoBBox.h:18
virtual const Double_t * GetOrigin() const
Definition TGeoBBox.h:77
virtual void InspectShape() const
Prints shape parameters.
Definition TGeoBBox.cxx:819
virtual Double_t GetDX() const
Definition TGeoBBox.h:74
virtual Double_t GetDZ() const
Definition TGeoBBox.h:76
virtual void Sizeof3D() const
virtual Double_t GetDY() const
Definition TGeoBBox.h:75
Double_t fOrigin[3]
Definition TGeoBBox.h:24
void SetBoxDimensions(Double_t dx, Double_t dy, Double_t dz, Double_t *origin=0)
Set parameters of the box.
Definition TGeoBBox.cxx:926
virtual void GetBoundingCylinder(Double_t *param) const
Fill vector param[4] with the bounding cylinder parameters.
Definition TGeoBBox.cxx:609
TGeoVolumeMulti * MakeVolumeMulti(const char *name, TGeoMedium *medium)
Make a TGeoVolumeMulti handling a list of volumes.
Geometrical transformation package.
Definition TGeoMatrix.h:41
Bool_t IsRotation() const
Definition TGeoMatrix.h:68
virtual void LocalToMaster(const Double_t *local, Double_t *master) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
Node containing an offset.
Definition TGeoNode.h:184
Parallelepiped class.
Definition TGeoPara.h:18
virtual void Sizeof3D() const
fill size of this 3-D object
Definition TGeoPara.cxx:661
virtual TGeoShape * GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix *mat) const
in case shape has some negative parameters, these has to be computed in order to fit the mother
Definition TGeoPara.cxx:526
Double_t fTxy
Definition TGeoPara.h:27
virtual Int_t GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const
Fills real parameters of a positioned box inside this. Returns 0 if successful.
Definition TGeoPara.cxx:452
virtual void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
Definition TGeoPara.cxx:689
Double_t fTyz
Definition TGeoPara.h:29
virtual void Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
Check the inside status for each of the points in the array.
Definition TGeoPara.cxx:671
Double_t fX
Definition TGeoPara.h:21
virtual Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0) const
compute distance from inside point to surface of the para Boundary safe algorithm.
Definition TGeoPara.cxx:244
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Definition TGeoPara.cxx:588
virtual void SetDimensions(Double_t *param)
Set dimensions starting from an array.
Definition TGeoPara.cxx:605
virtual void GetBoundingCylinder(Double_t *param) const
Fill vector param[4] with the bounding cylinder parameters.
Definition TGeoPara.cxx:444
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
Divide this parallelepiped shape belonging to volume "voldiv" into ndiv equal volumes called divname,...
Definition TGeoPara.cxx:371
TGeoPara()
Default constructor.
Definition TGeoPara.cxx:72
virtual void ComputeBBox()
compute bounding box
Definition TGeoPara.cxx:173
virtual void SetPoints(Double_t *points) const
Create PARA mesh points.
Definition TGeoPara.cxx:623
virtual Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const
computes the closest distance from given point to this shape, according to option.
Definition TGeoPara.cxx:563
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
Compute normal to closest surface from POINT.
Definition TGeoPara.cxx:185
Double_t fTheta
Definition TGeoPara.h:25
virtual ~TGeoPara()
destructor
Definition TGeoPara.cxx:157
Double_t fTxz
Definition TGeoPara.h:28
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
Definition TGeoPara.cxx:164
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const
Get range of shape for a given axis.
Definition TGeoPara.cxx:415
Double_t fAlpha
Definition TGeoPara.h:24
virtual void InspectShape() const
print shape parameters
Definition TGeoPara.cxx:546
virtual void DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
Definition TGeoPara.cxx:697
virtual void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
Compute the normal for an array o points so that norm.dot.dir is positive Input: Arrays of point coor...
Definition TGeoPara.cxx:681
virtual void Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
Compute safe distance from each of the points in the input array.
Definition TGeoPara.cxx:707
Double_t fZ
Definition TGeoPara.h:23
virtual Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0) const
compute distance from inside point to surface of the para
Definition TGeoPara.cxx:288
Double_t fY
Definition TGeoPara.h:22
Double_t fPhi
Definition TGeoPara.h:26
virtual Bool_t Contains(const Double_t *point) const
test if point is inside this sphere test Z range
Definition TGeoPara.cxx:229
Base finder class for patterns.
void SetDivIndex(Int_t index)
Base abstract class for all shapes.
Definition TGeoShape.h:26
static Double_t Big()
Definition TGeoShape.h:88
void SetShapeBit(UInt_t f, Bool_t set)
Equivalent of TObject::SetBit.
static Bool_t IsSameWithinTolerance(Double_t a, Double_t b)
Check if two numbers differ with less than a tolerance.
const char * GetPointerName() const
Provide a pointer name containing uid.
virtual const char * GetName() const
Get the shape name.
@ kGeoSavePrimitive
Definition TGeoShape.h:65
@ kGeoRunTimeShape
Definition TGeoShape.h:41
Bool_t TestShapeBit(UInt_t f) const
Definition TGeoShape.h:163
Volume families.
Definition TGeoVolume.h:255
void AddVolume(TGeoVolume *vol)
Add a volume with valid shape to the list of volumes.
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:49
void AddNodeOffset(TGeoVolume *vol, Int_t copy_no, Double_t offset=0, Option_t *option="")
Add a division node to the list of nodes.
TGeoMedium * GetMedium() const
Definition TGeoVolume.h:174
void SetFinder(TGeoPatternFinder *finder)
Definition TGeoVolume.h:232
Int_t GetNdaughters() const
Definition TGeoVolume.h:351
TObjArray * GetNodes()
Definition TGeoVolume.h:168
TObject * At(Int_t idx) const
Definition TObjArray.h:164
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:201
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:200
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:766
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:963
Basic string class.
Definition TString.h:136
const char * Data() const
Definition TString.h:369
Long64_t LocMin(Long64_t n, const T *a)
Return index of array with the minimum element.
Definition TMath.h:922
T1 Sign(T1 a, T2 b)
Definition TMathBase.h:161
Long64_t LocMax(Long64_t n, const T *a)
Return index of array with the maximum element.
Definition TMath.h:950
constexpr Double_t DegToRad()
Conversion from degree to radian:
Definition TMath.h:81
Double_t Sqrt(Double_t x)
Definition TMath.h:641
Short_t Min(Short_t a, Short_t b)
Definition TMathBase.h:176
Double_t Cos(Double_t)
Definition TMath.h:593
Double_t Sin(Double_t)
Definition TMath.h:589
Double_t Tan(Double_t)
Definition TMath.h:597
Short_t Abs(Short_t d)
Definition TMathBase.h:120