source: trunk/source/geometry/divisions/src/G4ParameterisationTubs.cc@ 1350

Last change on this file since 1350 was 1347, checked in by garnier, 15 years ago

geant4 tag 9.4

File size: 12.6 KB
RevLine 
[831]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//
[1347]27// $Id: G4ParameterisationTubs.cc,v 1.10 2010/11/10 09:16:13 gcosmo Exp $
28// GEANT4 tag $Name: geant4-09-04-ref-00 $
[831]29//
30// class G4ParameterisationTubs Implementation file
31//
32// 26.05.03 - P.Arce, Initial version
33// 08.04.04 - I.Hrivnacova, Implemented reflection
[1347]34// 21.04.10 - M.Asai, Added gaps
[831]35// --------------------------------------------------------------------
36
37#include "G4ParameterisationTubs.hh"
38
39#include <iomanip>
40#include "G4ThreeVector.hh"
41#include "G4RotationMatrix.hh"
42#include "G4VPhysicalVolume.hh"
43#include "G4LogicalVolume.hh"
44#include "G4ReflectedSolid.hh"
45#include "G4Tubs.hh"
46
47//--------------------------------------------------------------------------
48G4VParameterisationTubs::
49G4VParameterisationTubs( EAxis axis, G4int nDiv, G4double width,
50 G4double offset, G4VSolid* msolid,
51 DivisionType divType )
52 : G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid )
53{
54 G4Tubs* msol = (G4Tubs*)(msolid);
55 if (msolid->GetEntityType() == "G4ReflectedSolid")
56 {
57 //----- get constituent solid
58 G4VSolid* mConstituentSolid
59 = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
60 msol = (G4Tubs*)(mConstituentSolid);
61 fmotherSolid = msol;
62 fReflectedSolid = true;
63 }
64}
65
66//------------------------------------------------------------------------
67G4VParameterisationTubs::~G4VParameterisationTubs()
68{
69}
70
71//--------------------------------------------------------------------------
72G4ParameterisationTubsRho::
73G4ParameterisationTubsRho( EAxis axis, G4int nDiv,
74 G4double width, G4double offset,
75 G4VSolid* msolid, DivisionType divType )
76 : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
77{
78 CheckParametersValidity();
79 SetType( "DivisionTubsRho" );
80
81 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
82 if( divType == DivWIDTH )
83 {
84 fnDiv = CalculateNDiv( msol->GetOuterRadius() - msol->GetInnerRadius(),
85 width, offset );
86 }
87 else if( divType == DivNDIV )
88 {
89 fwidth = CalculateWidth( msol->GetOuterRadius() - msol->GetInnerRadius(),
90 nDiv, offset );
91 }
92
93#ifdef G4DIVDEBUG
94 if( verbose >= 1 )
95 {
96 G4cout << " G4ParameterisationTubsRho - no divisions " << fnDiv << " = "
97 << nDiv << G4endl
98 << " Offset " << foffset << " = " << offset << G4endl
99 << " Width " << fwidth << " = " << width << G4endl
100 << " DivType " << divType << G4endl;
101 }
102#endif
103}
104
105//--------------------------------------------------------------------------
106G4ParameterisationTubsRho::~G4ParameterisationTubsRho()
107{
108}
109
110//------------------------------------------------------------------------
111G4double G4ParameterisationTubsRho::GetMaxParameter() const
112{
113 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
114 return msol->GetOuterRadius() - msol->GetInnerRadius();
115}
116
117
118//--------------------------------------------------------------------------
119void
120G4ParameterisationTubsRho::
121ComputeTransformation(const G4int, G4VPhysicalVolume* physVol) const
122{
123 //----- translation
124 G4ThreeVector origin(0.,0.,0.);
125 //----- set translation
126 physVol->SetTranslation( origin );
127
128 //----- calculate rotation matrix: unit
129
130#ifdef G4DIVDEBUG
131 if( verbose >= 2 )
132 {
133 G4cout << " G4ParameterisationTubsRho " << G4endl
134 << " Offset: " << foffset/deg
135 << " - Width: " << fwidth/deg << G4endl;
136 }
137#endif
138
139 ChangeRotMatrix( physVol );
140
141#ifdef G4DIVDEBUG
142 if( verbose >= 2 )
143 {
144 G4cout << std::setprecision(8) << " G4ParameterisationTubsRho " << G4endl
145 << " Position: " << origin << " - Width: " << fwidth
146 << " - Axis " << faxis << G4endl;
147 }
148#endif
149}
150
151//--------------------------------------------------------------------------
152void
153G4ParameterisationTubsRho::
154ComputeDimensions( G4Tubs& tubs, const G4int copyNo,
155 const G4VPhysicalVolume* ) const
156{
157 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
158
[1347]159 G4double pRMin = msol->GetInnerRadius() + foffset + fwidth * copyNo + fhgap;
160 G4double pRMax = msol->GetInnerRadius() + foffset + fwidth * (copyNo+1) - fhgap;
[831]161 G4double pDz = msol->GetZHalfLength();
162 //- already rotated G4double pSR = foffset + copyNo*fwidth;
163 G4double pSPhi = msol->GetStartPhiAngle();
164 G4double pDPhi = msol->GetDeltaPhiAngle();;
165
166 tubs.SetInnerRadius( pRMin );
167 tubs.SetOuterRadius( pRMax );
168 tubs.SetZHalfLength( pDz );
169 tubs.SetStartPhiAngle( pSPhi );
170 tubs.SetDeltaPhiAngle( pDPhi );
171
172#ifdef G4DIVDEBUG
173 if( verbose >= 2 )
174 {
175 G4cout << " G4ParameterisationTubsRho::ComputeDimensions()" << G4endl
176 << " pRMin: " << pRMin << " - pRMax: " << pRMax << G4endl;
177 tubs.DumpInfo();
178 }
179#endif
180}
181
182//--------------------------------------------------------------------------
183G4ParameterisationTubsPhi::
184G4ParameterisationTubsPhi( EAxis axis, G4int nDiv,
185 G4double width, G4double offset,
186 G4VSolid* msolid, DivisionType divType )
187 : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
188{
189 CheckParametersValidity();
190 SetType( "DivisionTubsPhi" );
191
192 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
193 if( divType == DivWIDTH )
194 {
195 fnDiv = CalculateNDiv( msol->GetDeltaPhiAngle(), width, offset );
196 }
197 else if( divType == DivNDIV )
198 {
199 fwidth = CalculateWidth( msol->GetDeltaPhiAngle(), nDiv, offset );
200 }
201
202#ifdef G4DIVDEBUG
203 if( verbose >= 1 )
204 {
205 G4cout << " G4ParameterisationTubsPhi no divisions " << fnDiv << " = "
206 << nDiv << G4endl
207 << " Offset " << foffset << " = " << offset << G4endl
208 << " Width " << fwidth << " = " << width << G4endl;
209 }
210#endif
211}
212
213//--------------------------------------------------------------------------
214G4ParameterisationTubsPhi::~G4ParameterisationTubsPhi()
215{
216}
217
218//------------------------------------------------------------------------
219G4double G4ParameterisationTubsPhi::GetMaxParameter() const
220{
221 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
222 return msol->GetDeltaPhiAngle();
223}
224
225//--------------------------------------------------------------------------
226void
227G4ParameterisationTubsPhi::
228ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
229{
230 //----- translation
231 G4ThreeVector origin(0.,0.,0.);
232 //----- set translation
233 physVol->SetTranslation( origin );
234
235 //----- calculate rotation matrix (so that all volumes point to the centre)
236 G4double posi = foffset + copyNo*fwidth;
237
238#ifdef G4DIVDEBUG
239 if( verbose >= 2 )
240 {
241 G4cout << " G4ParameterisationTubsPhi - position: " << posi/deg << G4endl
242 << " copyNo: " << copyNo << " - foffset: " << foffset/deg
243 << " - fwidth: " << fwidth/deg << G4endl;
244 }
245#endif
246
247 ChangeRotMatrix( physVol, -posi );
248
249#ifdef G4DIVDEBUG
250 if( verbose >= 2 )
251 {
252 G4cout << std::setprecision(8) << " G4ParameterisationTubsPhi " << copyNo
253 << G4endl
254 << " Position: " << origin << " - Width: " << fwidth
255 << " - Axis: " << faxis << G4endl;
256 }
257#endif
258}
259
260//--------------------------------------------------------------------------
261void
262G4ParameterisationTubsPhi::
263ComputeDimensions( G4Tubs& tubs, const G4int,
264 const G4VPhysicalVolume* ) const
265{
266 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
267
268 G4double pRMin = msol->GetInnerRadius();
269 G4double pRMax = msol->GetOuterRadius();
270 G4double pDz = msol->GetZHalfLength();
271 //----- already rotated in 'ComputeTransformation'
[1347]272 G4double pSPhi = msol->GetStartPhiAngle() + fhgap;
273 G4double pDPhi = fwidth - 2.*fhgap;
[831]274
275 tubs.SetInnerRadius( pRMin );
276 tubs.SetOuterRadius( pRMax );
277 tubs.SetZHalfLength( pDz );
278 tubs.SetStartPhiAngle( pSPhi );
279 tubs.SetDeltaPhiAngle( pDPhi );
280
281#ifdef G4DIVDEBUG
282 if( verbose >= 2 )
283 {
284 G4cout << " G4ParameterisationTubsPhi::ComputeDimensions" << G4endl
285 << " pSPhi: " << pSPhi << " - pDPhi: " << pDPhi << G4endl;
286 tubs.DumpInfo();
287 }
288#endif
289}
290
291//--------------------------------------------------------------------------
292G4ParameterisationTubsZ::
293G4ParameterisationTubsZ( EAxis axis, G4int nDiv,
294 G4double width, G4double offset,
295 G4VSolid* msolid, DivisionType divType )
296 : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
297{
298 CheckParametersValidity();
299 SetType( "DivisionTubsZ" );
300
301 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
302 if( divType == DivWIDTH )
303 {
304 fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(), width, offset );
305 }
306 else if( divType == DivNDIV )
307 {
308 fwidth = CalculateWidth( 2*msol->GetZHalfLength(), nDiv, offset );
309 }
310
311#ifdef G4DIVDEBUG
312 if( verbose >= 1 )
313 {
314 G4cout << " G4ParameterisationTubsZ: # divisions " << fnDiv << " = "
315 << nDiv << G4endl
316 << " Offset " << foffset << " = " << offset << G4endl
317 << " Width " << fwidth << " = " << width << G4endl;
318 }
319#endif
320}
321
322//--------------------------------------------------------------------------
323G4ParameterisationTubsZ::~G4ParameterisationTubsZ()
324{
325}
326
327//------------------------------------------------------------------------
328G4double G4ParameterisationTubsZ::GetMaxParameter() const
329{
330 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
331 return 2*msol->GetZHalfLength();
332}
333
334//--------------------------------------------------------------------------
335void
336G4ParameterisationTubsZ::
337ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
338{
339 //----- set translation: along Z axis
340 G4Tubs* motherTubs = (G4Tubs*)(fmotherSolid);
341 G4double posi = - motherTubs->GetZHalfLength() + OffsetZ()
342 + fwidth/2 + copyNo*fwidth;
343 G4ThreeVector origin(0.,0.,posi);
344 physVol->SetTranslation( origin );
345
346 //----- calculate rotation matrix: unit
347
348#ifdef G4DIVDEBUG
349 if( verbose >= 2 )
350 {
351 G4cout << " G4ParameterisationTubsZ::ComputeTransformation()" << G4endl
352 << " Position: " << posi << " - copyNo: " << copyNo << G4endl
353 << " foffset " << foffset/deg << " - fwidth " << fwidth/deg
354 << G4endl;
355 }
356#endif
357
358 ChangeRotMatrix( physVol );
359
360#ifdef G4DIVDEBUG
361 if( verbose >= 2 )
362 {
363 G4cout << std::setprecision(8) << " G4ParameterisationTubsZ " << copyNo
364 << G4endl
365 << " Position: " << origin << " - Width: " << fwidth
366 << " - Axis: " << faxis << G4endl;
367 }
368#endif
369}
370
371//--------------------------------------------------------------------------
372void
373G4ParameterisationTubsZ::
374ComputeDimensions( G4Tubs& tubs, const G4int,
375 const G4VPhysicalVolume* ) const
376{
377 G4Tubs* msol = (G4Tubs*)(fmotherSolid);
378
379 G4double pRMin = msol->GetInnerRadius();
380 G4double pRMax = msol->GetOuterRadius();
381 // G4double pDz = msol->GetZHalfLength() / GetNoDiv();
[1347]382 G4double pDz = fwidth/2. - fhgap;
[831]383 G4double pSPhi = msol->GetStartPhiAngle();
384 G4double pDPhi = msol->GetDeltaPhiAngle();
385
386 tubs.SetInnerRadius( pRMin );
387 tubs.SetOuterRadius( pRMax );
388 tubs.SetZHalfLength( pDz );
389 tubs.SetStartPhiAngle( pSPhi );
390 tubs.SetDeltaPhiAngle( pDPhi );
391
392#ifdef G4DIVDEBUG
393 if( verbose >= 2 )
394 {
395 G4cout << " G4ParameterisationTubsZ::ComputeDimensions()" << G4endl
396 << " pDz: " << pDz << G4endl;
397 tubs.DumpInfo();
398 }
399#endif
400}
Note: See TracBrowser for help on using the repository browser.