Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TPARA.cxx
Go to the documentation of this file.
1// @(#)root/g3d:$Id$
2// Author: Nenad Buncic 19/09/95
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 "TPARA.h"
13#include "TNode.h"
14#include "TMath.h"
15
17
18/** \class TPARA
19\ingroup g3d
20A parallelepiped.
21
22\image html g3d_para.png
23It has 9 parameters:
24
25 - name: name of the shape
26 - title: shape's title
27 - material: (see TMaterial)
28 - dx: half-length in x
29 - dy: half-length in y
30 - dz: half-length in z
31 - alpha: angle formed by the y axis and by the plane joining the
32 centre of the faces parallel to the z-x plane at -DY and +DY
33 - theta: polar angle of the line joining the centres of the faces
34 at -DZ and +DZ in z
35 - phi: azimuthal angle of the line joining the centres of the
36 faces at -DZ and +DZ in z
37*/
38
39////////////////////////////////////////////////////////////////////////////////
40/// PARA shape default constructor
41
43{
44 fAlpha = 0.;
45 fTheta = 0.;
46 fPhi = 0.;
47}
48
49////////////////////////////////////////////////////////////////////////////////
50/// PARA shape normal constructor
51
52TPARA::TPARA(const char *name, const char *title, const char *material, Float_t dx, Float_t dy, Float_t dz,
53 Float_t alpha, Float_t theta, Float_t phi) : TBRIK(name, title,material, dx, dy, dz)
54{
55 fAlpha = alpha;
56 fTheta = theta;
57 fPhi = phi;
58}
59
60////////////////////////////////////////////////////////////////////////////////
61/// PARA shape default destructor
62
64{
65}
66
67////////////////////////////////////////////////////////////////////////////////
68/// Create PARA points
69
71{
72 if (!points) return;
73 Float_t dx, dy, dz, theta, phi, alpha;
74 const Float_t pi = Float_t (TMath::Pi());
75
76 alpha = fAlpha * pi/180.0;
77 theta = fTheta * pi/180.0;
78 phi = fPhi * pi/180.0;
79
80 dx = TBRIK::fDx;
81 dy = TBRIK::fDy;
82 dz = TBRIK::fDz;
83
84 // Parallelepiped change angles to tangents (by Pavel Nevski 12/04/99)
85 Double_t txy = TMath::Tan(alpha);
86 Double_t tth = TMath::Tan(theta);
87 Double_t txz = tth*TMath::Cos(phi);
88 Double_t tyz = tth*TMath::Sin(phi);
89
90 *points++ = -dz*txz-txy*dy-dx ; *points++ = -dy-dz*tyz ; *points++ = -dz;
91 *points++ = -dz*txz+txy*dy-dx ; *points++ = +dy-dz*tyz ; *points++ = -dz; //3
92 *points++ = -dz*txz+txy*dy+dx ; *points++ = +dy-dz*tyz ; *points++ = -dz;
93 *points++ = -dz*txz-txy*dy+dx ; *points++ = -dy-dz*tyz ; *points++ = -dz;//1
94 *points++ = +dz*txz-txy*dy-dx ; *points++ = -dy+dz*tyz ; *points++ = +dz;
95 *points++ = +dz*txz+txy*dy-dx ; *points++ = +dy+dz*tyz ; *points++ = +dz;//7
96 *points++ = +dz*txz+txy*dy+dx ; *points++ = +dy+dz*tyz ; *points++ = +dz;
97 *points++ = +dz*txz-txy*dy+dx ; *points++ = -dy+dz*tyz ; *points++ = +dz;//5
98}
float Float_t
Definition RtypesCore.h:57
#define ClassImp(name)
Definition Rtypes.h:382
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
char name[80]
Definition TGX11.cxx:110
A box with faces perpendicular to the axes.
Definition TBRIK.h:26
Float_t fDz
Definition TBRIK.h:31
Float_t fDy
Definition TBRIK.h:30
Float_t fDx
Definition TBRIK.h:29
A parallelepiped.
Definition TPARA.h:30
Float_t fAlpha
Definition TPARA.h:32
Float_t fTheta
Definition TPARA.h:33
~TPARA() override
PARA shape default destructor.
Definition TPARA.cxx:63
void SetPoints(Double_t *points) const override
Create PARA points.
Definition TPARA.cxx:70
TPARA()
PARA shape default constructor.
Definition TPARA.cxx:42
Float_t fPhi
Definition TPARA.h:34
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:594
constexpr Double_t Pi()
Definition TMath.h:37
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:588
Double_t Tan(Double_t)
Returns the tangent of an angle of x radians.
Definition TMath.h:600