Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParameterisationTubs.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4ParameterisationTubs[Rho/Phi/Z] implementation
27//
28// 26.05.03 - P.Arce, Initial version
29// 08.04.04 - I.Hrivnacova, Implemented reflection
30// 21.04.10 - M.Asai, Added gaps
31// --------------------------------------------------------------------
32
34
35#include <iomanip>
36#include "G4ThreeVector.hh"
37#include "G4RotationMatrix.hh"
38#include "G4VPhysicalVolume.hh"
39#include "G4LogicalVolume.hh"
40#include "G4ReflectedSolid.hh"
41#include "G4Tubs.hh"
42
43//--------------------------------------------------------------------------
46 G4double offset, G4VSolid* msolid,
47 DivisionType divType )
48 : G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid )
49{
50 G4Tubs* msol = (G4Tubs*)(msolid);
51 if (msolid->GetEntityType() == "G4ReflectedSolid")
52 {
53 //----- get constituent solid
54 G4VSolid* mConstituentSolid
55 = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
56 msol = (G4Tubs*)(mConstituentSolid);
57 fmotherSolid = msol;
58 fReflectedSolid = true;
59 }
60}
61
62//------------------------------------------------------------------------
64{
65}
66
67//--------------------------------------------------------------------------
70 G4double width, G4double offset,
71 G4VSolid* msolid, DivisionType divType )
72 : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
73{
75 SetType( "DivisionTubsRho" );
76
77 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
78 if( divType == DivWIDTH )
79 {
81 width, offset );
82 }
83 else if( divType == DivNDIV )
84 {
86 nDiv, offset );
87 }
88
89#ifdef G4DIVDEBUG
90 if( verbose >= 1 )
91 {
92 G4cout << " G4ParameterisationTubsRho - no divisions " << fnDiv << " = "
93 << nDiv << G4endl
94 << " Offset " << foffset << " = " << offset << G4endl
95 << " Width " << fwidth << " = " << width << G4endl
96 << " DivType " << divType << G4endl;
97 }
98#endif
99}
100
101//--------------------------------------------------------------------------
103{
104}
105
106//------------------------------------------------------------------------
108{
109 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
110 return msol->GetOuterRadius() - msol->GetInnerRadius();
111}
112
113
114//--------------------------------------------------------------------------
115void
117ComputeTransformation(const G4int, G4VPhysicalVolume* physVol) const
118{
119 //----- translation
120 G4ThreeVector origin(0.,0.,0.);
121 //----- set translation
122 physVol->SetTranslation( origin );
123
124 //----- calculate rotation matrix: unit
125
126#ifdef G4DIVDEBUG
127 if( verbose >= 2 )
128 {
129 G4cout << " G4ParameterisationTubsRho " << G4endl
130 << " Offset: " << foffset/CLHEP::deg
131 << " - Width: " << fwidth/CLHEP::deg << G4endl;
132 }
133#endif
134
135 ChangeRotMatrix( physVol );
136
137#ifdef G4DIVDEBUG
138 if( verbose >= 2 )
139 {
140 G4cout << std::setprecision(8) << " G4ParameterisationTubsRho " << G4endl
141 << " Position: " << origin << " - Width: " << fwidth
142 << " - Axis " << faxis << G4endl;
143 }
144#endif
145}
146
147//--------------------------------------------------------------------------
148void
150ComputeDimensions( G4Tubs& tubs, const G4int copyNo,
151 const G4VPhysicalVolume* ) const
152{
153 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
154
155 G4double pRMin = msol->GetInnerRadius() + foffset + fwidth*copyNo + fhgap;
156 G4double pRMax = msol->GetInnerRadius() + foffset + fwidth*(copyNo+1) - fhgap;
157 G4double pDz = msol->GetZHalfLength();
158 //- already rotated G4double pSR = foffset + copyNo*fwidth;
159 G4double pSPhi = msol->GetStartPhiAngle();
160 G4double pDPhi = msol->GetDeltaPhiAngle();;
161
162 tubs.SetInnerRadius( pRMin );
163 tubs.SetOuterRadius( pRMax );
164 tubs.SetZHalfLength( pDz );
165 tubs.SetStartPhiAngle( pSPhi, false );
166 tubs.SetDeltaPhiAngle( pDPhi );
167
168#ifdef G4DIVDEBUG
169 if( verbose >= 2 )
170 {
171 G4cout << " G4ParameterisationTubsRho::ComputeDimensions()" << G4endl
172 << " pRMin: " << pRMin << " - pRMax: " << pRMax << G4endl;
173 tubs.DumpInfo();
174 }
175#endif
176}
177
178//--------------------------------------------------------------------------
181 G4double width, G4double offset,
182 G4VSolid* msolid, DivisionType divType )
183 : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
184{
186 SetType( "DivisionTubsPhi" );
187
188 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
189 if( divType == DivWIDTH )
190 {
191 fnDiv = CalculateNDiv( msol->GetDeltaPhiAngle(), width, offset );
192 }
193 else if( divType == DivNDIV )
194 {
195 fwidth = CalculateWidth( msol->GetDeltaPhiAngle(), nDiv, offset );
196 }
197
198#ifdef G4DIVDEBUG
199 if( verbose >= 1 )
200 {
201 G4cout << " G4ParameterisationTubsPhi no divisions " << fnDiv << " = "
202 << nDiv << G4endl
203 << " Offset " << foffset << " = " << offset << G4endl
204 << " Width " << fwidth << " = " << width << G4endl;
205 }
206#endif
207}
208
209//--------------------------------------------------------------------------
211{
212}
213
214//------------------------------------------------------------------------
216{
217 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
218 return msol->GetDeltaPhiAngle();
219}
220
221//--------------------------------------------------------------------------
222void
224ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
225{
226 //----- translation
227 G4ThreeVector origin(0.,0.,0.);
228 //----- set translation
229 physVol->SetTranslation( origin );
230
231 //----- calculate rotation matrix (so that all volumes point to the centre)
232 G4double posi = foffset + copyNo*fwidth;
233
234#ifdef G4DIVDEBUG
235 if( verbose >= 2 )
236 {
237 G4cout << " G4ParameterisationTubsPhi - position: " << posi/CLHEP::deg << G4endl
238 << " copyNo: " << copyNo << " - foffset: " << foffset/CLHEP::deg
239 << " - fwidth: " << fwidth/CLHEP::deg << G4endl;
240 }
241#endif
242
243 ChangeRotMatrix( physVol, -posi );
244
245#ifdef G4DIVDEBUG
246 if( verbose >= 2 )
247 {
248 G4cout << std::setprecision(8) << " G4ParameterisationTubsPhi " << copyNo
249 << G4endl
250 << " Position: " << origin << " - Width: " << fwidth
251 << " - Axis: " << faxis << G4endl;
252 }
253#endif
254}
255
256//--------------------------------------------------------------------------
257void
259ComputeDimensions( G4Tubs& tubs, const G4int,
260 const G4VPhysicalVolume* ) const
261{
262 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
263
264 G4double pRMin = msol->GetInnerRadius();
265 G4double pRMax = msol->GetOuterRadius();
266 G4double pDz = msol->GetZHalfLength();
267 //----- already rotated in 'ComputeTransformation'
268 G4double pSPhi = msol->GetStartPhiAngle() + fhgap;
269 G4double pDPhi = fwidth - 2.*fhgap;
270
271 tubs.SetInnerRadius( pRMin );
272 tubs.SetOuterRadius( pRMax );
273 tubs.SetZHalfLength( pDz );
274 tubs.SetStartPhiAngle( pSPhi, false );
275 tubs.SetDeltaPhiAngle( pDPhi );
276
277#ifdef G4DIVDEBUG
278 if( verbose >= 2 )
279 {
280 G4cout << " G4ParameterisationTubsPhi::ComputeDimensions" << G4endl
281 << " pSPhi: " << pSPhi << " - pDPhi: " << pDPhi << G4endl;
282 tubs.DumpInfo();
283 }
284#endif
285}
286
287//--------------------------------------------------------------------------
290 G4double width, G4double offset,
291 G4VSolid* msolid, DivisionType divType )
292 : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
293{
295 SetType( "DivisionTubsZ" );
296
297 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
298 if( divType == DivWIDTH )
299 {
300 fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(), width, offset );
301 }
302 else if( divType == DivNDIV )
303 {
304 fwidth = CalculateWidth( 2*msol->GetZHalfLength(), nDiv, offset );
305 }
306
307#ifdef G4DIVDEBUG
308 if( verbose >= 1 )
309 {
310 G4cout << " G4ParameterisationTubsZ: # divisions " << fnDiv << " = "
311 << nDiv << G4endl
312 << " Offset " << foffset << " = " << offset << G4endl
313 << " Width " << fwidth << " = " << width << G4endl;
314 }
315#endif
316}
317
318//--------------------------------------------------------------------------
320{
321}
322
323//------------------------------------------------------------------------
325{
326 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
327 return 2*msol->GetZHalfLength();
328}
329
330//--------------------------------------------------------------------------
331void
333ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
334{
335 //----- set translation: along Z axis
336 G4Tubs* motherTubs = (G4Tubs*)(fmotherSolid);
337 G4double posi = - motherTubs->GetZHalfLength() + OffsetZ()
338 + fwidth/2 + copyNo*fwidth;
339 G4ThreeVector origin(0.,0.,posi);
340 physVol->SetTranslation( origin );
341
342 //----- calculate rotation matrix: unit
343
344#ifdef G4DIVDEBUG
345 if( verbose >= 2 )
346 {
347 G4cout << " G4ParameterisationTubsZ::ComputeTransformation()" << G4endl
348 << " Position: " << posi << " - copyNo: " << copyNo << G4endl
349 << " foffset " << foffset/CLHEP::deg << " - fwidth " << fwidth/CLHEP::deg
350 << G4endl;
351 }
352#endif
353
354 ChangeRotMatrix( physVol );
355
356#ifdef G4DIVDEBUG
357 if( verbose >= 2 )
358 {
359 G4cout << std::setprecision(8) << " G4ParameterisationTubsZ " << copyNo
360 << G4endl
361 << " Position: " << origin << " - Width: " << fwidth
362 << " - Axis: " << faxis << G4endl;
363 }
364#endif
365}
366
367//--------------------------------------------------------------------------
368void
370ComputeDimensions( G4Tubs& tubs, const G4int,
371 const G4VPhysicalVolume* ) const
372{
373 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
374
375 G4double pRMin = msol->GetInnerRadius();
376 G4double pRMax = msol->GetOuterRadius();
377 // G4double pDz = msol->GetZHalfLength() / GetNoDiv();
378 G4double pDz = fwidth/2. - fhgap;
379 G4double pSPhi = msol->GetStartPhiAngle();
380 G4double pDPhi = msol->GetDeltaPhiAngle();
381
382 tubs.SetInnerRadius( pRMin );
383 tubs.SetOuterRadius( pRMax );
384 tubs.SetZHalfLength( pDz );
385 tubs.SetStartPhiAngle( pSPhi, false );
386 tubs.SetDeltaPhiAngle( pDPhi );
387
388#ifdef G4DIVDEBUG
389 if( verbose >= 2 )
390 {
391 G4cout << " G4ParameterisationTubsZ::ComputeDimensions()" << G4endl
392 << " pDz: " << pDz << G4endl;
393 tubs.DumpInfo();
394 }
395#endif
396}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParameterisationTubsPhi(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
void ComputeDimensions(G4Tubs &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
void ComputeDimensions(G4Tubs &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
G4ParameterisationTubsRho(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
G4ParameterisationTubsZ(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
void ComputeDimensions(G4Tubs &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const
Definition: G4Tubs.hh:75
G4double GetZHalfLength() const
void SetDeltaPhiAngle(G4double newDPhi)
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
G4double GetInnerRadius() const
G4double GetOuterRadius() const
G4double GetStartPhiAngle() const
void SetInnerRadius(G4double newRMin)
void SetOuterRadius(G4double newRMax)
G4double GetDeltaPhiAngle() const
void SetZHalfLength(G4double newDz)
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.0) const
G4VParameterisationTubs(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
void SetTranslation(const G4ThreeVector &v)
void DumpInfo() const
virtual G4GeometryType GetEntityType() const =0
EAxis
Definition: geomdefs.hh:54