source: trunk/source/geometry/divisions/src/G4ParameterisationCons.cc@ 1279

Last change on this file since 1279 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 14.4 KB
Line 
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//
27// $Id: G4ParameterisationCons.cc,v 1.9 2006/06/29 18:18:38 gunter Exp $
28// GEANT4 tag $Name: geant4-09-03 $
29//
30// class G4ParameterisationCons Implementation file
31//
32// 26.05.03 - P.Arce, Initial version
33// 08.04.04 - I.Hrivnacova, Implemented reflection
34// --------------------------------------------------------------------
35
36#include "G4ParameterisationCons.hh"
37
38#include <iomanip>
39#include "G4ThreeVector.hh"
40#include "G4RotationMatrix.hh"
41#include "G4VPhysicalVolume.hh"
42#include "G4LogicalVolume.hh"
43#include "G4ReflectedSolid.hh"
44#include "G4Cons.hh"
45
46//--------------------------------------------------------------------------
47G4VParameterisationCons::
48G4VParameterisationCons( EAxis axis, G4int nDiv, G4double width,
49 G4double offset, G4VSolid* msolid,
50 DivisionType divType )
51 : G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid )
52{
53 G4Cons* msol = (G4Cons*)(msolid);
54 if (msolid->GetEntityType() == "G4ReflectedSolid")
55 {
56 // Get constituent solid
57 G4VSolid* mConstituentSolid
58 = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
59 msol = (G4Cons*)(mConstituentSolid);
60
61 // Create a new solid with inversed parameters
62 G4Cons* newSolid
63 = new G4Cons(msol->GetName(),
64 msol->GetInnerRadiusPlusZ(), msol->GetOuterRadiusPlusZ(),
65 msol->GetInnerRadiusMinusZ(), msol->GetOuterRadiusMinusZ(),
66 msol->GetZHalfLength(),
67 msol->GetStartPhiAngle(), msol->GetDeltaPhiAngle());
68 msol = newSolid;
69 fmotherSolid = newSolid;
70 fReflectedSolid = true;
71 fDeleteSolid = true;
72 }
73}
74
75//------------------------------------------------------------------------
76G4VParameterisationCons::~G4VParameterisationCons()
77{
78}
79
80//--------------------------------------------------------------------------
81G4ParameterisationConsRho::
82G4ParameterisationConsRho( EAxis axis, G4int nDiv,
83 G4double width, G4double offset,
84 G4VSolid* msolid, DivisionType divType )
85 : G4VParameterisationCons( axis, nDiv, width, offset, msolid, divType )
86{
87 CheckParametersValidity();
88 SetType( "DivisionConsRho" );
89
90 G4Cons* msol = (G4Cons*)(fmotherSolid);
91 if( msol->GetInnerRadiusPlusZ() == 0. )
92 {
93 G4cerr << "WARNING - G4ParameterisationConsRho, OuterRadiusMinusZ = 0. "
94 << G4endl
95 << " Width is calculated as that of OuterRadiusMinusZ !"
96 << G4endl;
97 }
98
99 if( divType == DivWIDTH )
100 {
101 fnDiv = CalculateNDiv( msol->GetOuterRadiusMinusZ()
102 - msol->GetInnerRadiusMinusZ(), width, offset );
103 }
104 else if( divType == DivNDIV )
105 {
106 G4Cons* msol = (G4Cons*)(msolid);
107 fwidth = CalculateWidth( msol->GetOuterRadiusMinusZ()
108 - msol->GetInnerRadiusMinusZ(), nDiv, offset );
109 }
110
111#ifdef G4DIVDEBUG
112 if( verbose >= 1 )
113 {
114 G4cout << " G4ParameterisationConsRho - no divisions " << fnDiv << " = "
115 << nDiv << G4endl
116 << " Offset " << foffset << " = " << offset
117 << " - Width " << fwidth << " = " << width << G4endl;
118 }
119#endif
120}
121
122//--------------------------------------------------------------------------
123G4ParameterisationConsRho::~G4ParameterisationConsRho()
124{
125}
126
127//------------------------------------------------------------------------
128G4double G4ParameterisationConsRho::GetMaxParameter() const
129{
130 G4Cons* msol = (G4Cons*)(fmotherSolid);
131 return msol->GetOuterRadiusMinusZ() - msol->GetInnerRadiusMinusZ();
132}
133
134//--------------------------------------------------------------------------
135void
136G4ParameterisationConsRho::
137ComputeTransformation( const G4int, G4VPhysicalVolume *physVol ) const
138{
139 //----- translation
140 G4ThreeVector origin(0.,0.,0.);
141 //----- set translation
142 physVol->SetTranslation( origin );
143
144 //----- calculate rotation matrix: unit
145
146#ifdef G4DIVDEBUG
147 if( verbose >= 2 )
148 {
149 G4cout << " G4ParameterisationConsRho " << G4endl
150 << " Offset: " << foffset
151 << " - Width: " << fwidth << G4endl;
152 }
153#endif
154
155 ChangeRotMatrix( physVol );
156
157#ifdef G4DIVDEBUG
158 if( verbose >= 2 )
159 {
160 G4cout << std::setprecision(8) << " G4ParameterisationConsRho" << G4endl
161 << " Position: " << origin << " - Width: " << fwidth
162 << " - Axis: " << faxis << G4endl;
163 }
164#endif
165}
166
167//--------------------------------------------------------------------------
168void
169G4ParameterisationConsRho::
170ComputeDimensions( G4Cons& cons, const G4int copyNo,
171 const G4VPhysicalVolume* ) const
172{
173 G4Cons* msol = (G4Cons*)(fmotherSolid);
174
175 G4double pRMin1 = msol->GetInnerRadiusMinusZ() + foffset + fwidth * copyNo;
176 G4double pRMax1 = msol->GetInnerRadiusMinusZ() + foffset + fwidth * (copyNo+1);
177
178 //width at Z Plus
179 //- G4double fwidthPlus =
180 // fwidth * ( msol->GetOuterRadiusPlusZ()/ msol->GetInnerRadiusPlusZ())
181 //- / ( msol->GetOuterRadiusMinusZ() - msol->GetInnerRadiusMinusZ());
182 G4double fwidthPlus = CalculateWidth( msol->GetOuterRadiusPlusZ()
183 - msol->GetInnerRadiusPlusZ(), fnDiv, foffset );
184 G4double pRMin2 = msol->GetInnerRadiusPlusZ()
185 + foffset + fwidthPlus * copyNo;
186 G4double pRMax2 = msol->GetInnerRadiusPlusZ()
187 + foffset + fwidthPlus * (copyNo+1);
188 G4double pDz = msol->GetZHalfLength();
189
190 //- already rotated double pSR = foffset + copyNo*fwidth;
191 G4double pSPhi = msol->GetStartPhiAngle();
192 G4double pDPhi = msol->GetDeltaPhiAngle();;
193
194 cons.SetInnerRadiusMinusZ( pRMin1 );
195 cons.SetOuterRadiusMinusZ( pRMax1 );
196 cons.SetInnerRadiusPlusZ( pRMin2 );
197 cons.SetOuterRadiusPlusZ( pRMax2 );
198 cons.SetZHalfLength( pDz );
199 cons.SetStartPhiAngle( pSPhi );
200 cons.SetDeltaPhiAngle( pDPhi );
201
202#ifdef G4DIVDEBUG
203 if( verbose >= 2 )
204 {
205 G4cout << " G4ParameterisationConsRho::ComputeDimensions()" << G4endl
206 << " pRMin: " << pRMin1 << " - pRMax: " << pRMax1 << G4endl;
207 if( verbose >= 4 ) cons.DumpInfo();
208 }
209#endif
210}
211
212//--------------------------------------------------------------------------
213G4ParameterisationConsPhi::
214G4ParameterisationConsPhi( EAxis axis, G4int nDiv,
215 G4double width, G4double offset,
216 G4VSolid* msolid, DivisionType divType )
217 : G4VParameterisationCons( axis, nDiv, width, offset, msolid, divType )
218{
219 CheckParametersValidity();
220 SetType( "DivisionConsPhi" );
221
222 G4Cons* msol = (G4Cons*)(fmotherSolid);
223 if( divType == DivWIDTH )
224 {
225 fnDiv = CalculateNDiv( msol->GetDeltaPhiAngle(), width, offset );
226 }
227 else if( divType == DivNDIV )
228 {
229 fwidth = CalculateWidth( msol->GetDeltaPhiAngle(), nDiv, offset );
230 }
231
232#ifdef G4DIVDEBUG
233 if( verbose >= 1 )
234 {
235 G4cout << " G4ParameterisationConsPhi no divisions " << fnDiv << " = "
236 << nDiv << G4endl
237 << " Offset " << foffset << " = " << offset << G4endl
238 << " Width " << fwidth << " = " << width << G4endl;
239 }
240#endif
241}
242
243//--------------------------------------------------------------------------
244G4ParameterisationConsPhi::~G4ParameterisationConsPhi()
245{
246}
247
248//------------------------------------------------------------------------
249G4double G4ParameterisationConsPhi::GetMaxParameter() const
250{
251 G4Cons* msol = (G4Cons*)(fmotherSolid);
252 return msol->GetDeltaPhiAngle();
253}
254
255//--------------------------------------------------------------------------
256void
257G4ParameterisationConsPhi::
258ComputeTransformation( const G4int copyNo, G4VPhysicalVolume *physVol ) const
259{
260 //----- translation
261 G4ThreeVector origin(0.,0.,0.);
262 //----- set translation
263 physVol->SetTranslation( origin );
264
265 //----- calculate rotation matrix (so that all volumes point to the centre)
266 G4double posi = foffset + copyNo*fwidth;
267
268#ifdef G4DIVDEBUG
269 if( verbose >= 2 )
270 {
271 G4cout << " G4ParameterisationConsPhi - position: " << posi/deg << G4endl
272 << " Origin: " << origin << " copyNo: " << copyNo
273 << " - foffset: " << foffset/deg
274 << " - fwidth: " << fwidth/deg << G4endl
275 << " - Axis: " << faxis << G4endl;
276 }
277#endif
278
279 ChangeRotMatrix( physVol, -posi );
280}
281
282//--------------------------------------------------------------------------
283void
284G4ParameterisationConsPhi::
285ComputeDimensions( G4Cons& cons, const G4int,
286 const G4VPhysicalVolume* ) const
287{
288 G4Cons* msol = (G4Cons*)(fmotherSolid);
289
290 G4double pRMin1 = msol->GetInnerRadiusMinusZ();
291 G4double pRMax1 = msol->GetOuterRadiusMinusZ();
292 G4double pRMin2 = msol->GetInnerRadiusPlusZ();
293 G4double pRMax2 = msol->GetOuterRadiusPlusZ();
294 G4double pDz = msol->GetZHalfLength();
295
296 //- already rotated double pSPhi = foffset + copyNo*fwidth;
297 G4double pSPhi = foffset + msol->GetStartPhiAngle();
298 G4double pDPhi = fwidth;
299
300 cons.SetInnerRadiusMinusZ( pRMin1 );
301 cons.SetOuterRadiusMinusZ( pRMax1 );
302 cons.SetInnerRadiusPlusZ( pRMin2 );
303 cons.SetOuterRadiusPlusZ( pRMax2 );
304 cons.SetZHalfLength( pDz );
305 cons.SetStartPhiAngle( pSPhi );
306 cons.SetDeltaPhiAngle( pDPhi );
307
308#ifdef G4DIVDEBUG
309 if( verbose >= 2 )
310 {
311 G4cout << " G4ParameterisationConsPhi::ComputeDimensions" << G4endl
312 << " pSPhi: " << pSPhi << " - pDPhi: " << pDPhi << G4endl;
313 if( verbose >= 4 ) cons.DumpInfo();
314 }
315#endif
316}
317
318//--------------------------------------------------------------------------
319G4ParameterisationConsZ::
320G4ParameterisationConsZ( EAxis axis, G4int nDiv,
321 G4double width, G4double offset,
322 G4VSolid* msolid, DivisionType divType )
323 : G4VParameterisationCons( axis, nDiv, width, offset, msolid, divType )
324{
325 CheckParametersValidity();
326 SetType( "DivisionConsZ" );
327
328 G4Cons* msol = (G4Cons*)(fmotherSolid);
329 if( divType == DivWIDTH )
330 {
331 fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(), width, offset );
332 }
333 else if( divType == DivNDIV )
334 {
335 fwidth = CalculateWidth( 2*msol->GetZHalfLength(), nDiv, offset );
336 }
337
338#ifdef G4DIVDEBUG
339 if( verbose >= 1 )
340 {
341 G4cout << " G4ParameterisationConsZ: # divisions " << fnDiv << " = "
342 << nDiv << G4endl
343 << " Offset " << foffset << " = " << offset << G4endl
344 << " Width " << fwidth << " = " << width << G4endl
345 << " - Axis: " << faxis << G4endl;
346 }
347#endif
348}
349
350//--------------------------------------------------------------------------
351G4ParameterisationConsZ::~G4ParameterisationConsZ()
352{
353}
354
355//------------------------------------------------------------------------
356G4double G4ParameterisationConsZ::GetMaxParameter() const
357{
358 G4Cons* msol = (G4Cons*)(fmotherSolid);
359 return 2*msol->GetZHalfLength();
360}
361
362//--------------------------------------------------------------------------
363void
364G4ParameterisationConsZ::
365ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol ) const
366{
367 //----- set translation: along Z axis
368 G4Cons* motherCons = (G4Cons*)(GetMotherSolid());
369 G4double posi = - motherCons->GetZHalfLength() + OffsetZ()
370 + fwidth/2 + copyNo*fwidth;
371 G4ThreeVector origin(0.,0.,posi);
372 physVol->SetTranslation( origin );
373
374 //----- calculate rotation matrix: unit
375
376#ifdef G4DIVDEBUG
377 if( verbose >= 2 )
378 {
379 G4cout << " G4ParameterisationConsZ::ComputeTransformation()" << G4endl
380 << " Origin: " << origin << " - copyNo: " << copyNo << G4endl
381 << " foffset: " << foffset << " - fwidth: " << fwidth
382 << G4endl;
383 }
384#endif
385
386 ChangeRotMatrix( physVol );
387}
388
389
390//--------------------------------------------------------------------------
391void
392G4ParameterisationConsZ::
393ComputeDimensions( G4Cons& cons, const G4int copyNo,
394 const G4VPhysicalVolume* ) const
395{
396 G4Cons* msol = (G4Cons*)(fmotherSolid);
397
398 G4double mHalfLength = msol->GetZHalfLength();
399 G4double aRInner = (msol->GetInnerRadiusPlusZ()
400 - msol->GetInnerRadiusMinusZ()) / (2*mHalfLength);
401 G4double bRInner = (msol->GetInnerRadiusPlusZ()
402 + msol->GetInnerRadiusMinusZ()) / 2;
403 G4double aROuter = (msol->GetOuterRadiusPlusZ()
404 - msol->GetOuterRadiusMinusZ()) / (2*mHalfLength);
405 G4double bROuter = (msol->GetOuterRadiusPlusZ()
406 + msol->GetOuterRadiusMinusZ()) / 2;
407 G4double xMinusZ = -mHalfLength + OffsetZ() + fwidth*copyNo;
408 G4double xPlusZ = -mHalfLength + OffsetZ() + fwidth*(copyNo+1);
409 cons.SetInnerRadiusMinusZ( aRInner * xMinusZ + bRInner );
410 cons.SetOuterRadiusMinusZ( aROuter * xMinusZ + bROuter );
411 cons.SetInnerRadiusPlusZ( aRInner * xPlusZ + bRInner );
412 cons.SetOuterRadiusPlusZ( aROuter * xPlusZ + bROuter );
413
414 G4double pDz = fwidth / 2.;
415 G4double pSPhi = msol->GetStartPhiAngle();
416 G4double pDPhi = msol->GetDeltaPhiAngle();
417
418 cons.SetZHalfLength( pDz );
419 cons.SetStartPhiAngle( pSPhi );
420 cons.SetDeltaPhiAngle( pDPhi );
421
422#ifdef G4DIVDEBUG
423 if( verbose >= 2 )
424 {
425 G4cout << " G4ParameterisationConsZ::ComputeDimensions()" << G4endl
426 << " pDz: " << pDz << G4endl;
427 if( verbose >= 4 ) cons.DumpInfo();
428 }
429#endif
430
431}
Note: See TracBrowser for help on using the repository browser.