Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParameterisationCons.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//
27// $Id$
28//
29// class G4ParameterisationCons Implementation file
30//
31// 26.05.03 - P.Arce, Initial version
32// 08.04.04 - I.Hrivnacova, Implemented reflection
33// 21.04.10 - M.Asai, Added gaps
34// --------------------------------------------------------------------
35
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//--------------------------------------------------------------------------
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(),
66 msol->GetZHalfLength(),
67 msol->GetStartPhiAngle(), msol->GetDeltaPhiAngle());
68 msol = newSolid;
69 fmotherSolid = newSolid;
70 fReflectedSolid = true;
71 fDeleteSolid = true;
72 }
73}
74
75//------------------------------------------------------------------------
77{
78}
79
80//--------------------------------------------------------------------------
83 G4double width, G4double offset,
84 G4VSolid* msolid, DivisionType divType )
85 : G4VParameterisationCons( axis, nDiv, width, offset, msolid, divType )
86{
88 SetType( "DivisionConsRho" );
89
90 G4Cons* msol = (G4Cons*)(fmotherSolid);
91 if( msol->GetInnerRadiusPlusZ() == 0. )
92 {
93 std::ostringstream message;
94 message << "OuterRadiusMinusZ = 0" << G4endl
95 << "Width is calculated as that of OuterRadiusMinusZ !";
96 G4Exception("G4ParameterisationConsRho::G4ParameterisationConsRho()",
97 "GeomDiv1001", JustWarning, message);
98 }
99
100 if( divType == DivWIDTH )
101 {
103 - msol->GetInnerRadiusMinusZ(), width, offset );
104 }
105 else if( divType == DivNDIV )
106 {
107 G4Cons* mconsol = (G4Cons*)(msolid);
109 - mconsol->GetInnerRadiusMinusZ(), nDiv, offset );
110 }
111
112#ifdef G4DIVDEBUG
113 if( verbose >= 1 )
114 {
115 G4cout << " G4ParameterisationConsRho - no divisions " << fnDiv << " = "
116 << nDiv << G4endl
117 << " Offset " << foffset << " = " << offset
118 << " - Width " << fwidth << " = " << width << G4endl;
119 }
120#endif
121}
122
123//--------------------------------------------------------------------------
125{
126}
127
128//------------------------------------------------------------------------
130{
131 G4Cons* msol = (G4Cons*)(fmotherSolid);
132 return msol->GetOuterRadiusMinusZ() - msol->GetInnerRadiusMinusZ();
133}
134
135//--------------------------------------------------------------------------
136void
138ComputeTransformation( const G4int, G4VPhysicalVolume *physVol ) const
139{
140 //----- translation
141 G4ThreeVector origin(0.,0.,0.);
142 //----- set translation
143 physVol->SetTranslation( origin );
144
145 //----- calculate rotation matrix: unit
146
147#ifdef G4DIVDEBUG
148 if( verbose >= 2 )
149 {
150 G4cout << " G4ParameterisationConsRho " << G4endl
151 << " Offset: " << foffset
152 << " - Width: " << fwidth << G4endl;
153 }
154#endif
155
156 ChangeRotMatrix( physVol );
157
158#ifdef G4DIVDEBUG
159 if( verbose >= 2 )
160 {
161 G4cout << std::setprecision(8) << " G4ParameterisationConsRho" << G4endl
162 << " Position: " << origin << " - Width: " << fwidth
163 << " - Axis: " << faxis << G4endl;
164 }
165#endif
166}
167
168//--------------------------------------------------------------------------
169void
171ComputeDimensions( G4Cons& cons, const G4int copyNo,
172 const G4VPhysicalVolume* ) const
173{
174 G4Cons* msol = (G4Cons*)(fmotherSolid);
175
176 G4double pRMin1 = msol->GetInnerRadiusMinusZ() + foffset + fwidth * copyNo;
177 G4double pRMax1 = msol->GetInnerRadiusMinusZ() + foffset + fwidth * (copyNo+1);
178
179 //width at Z Plus
180 //- G4double fwidthPlus =
181 // fwidth * ( msol->GetOuterRadiusPlusZ()/ msol->GetInnerRadiusPlusZ())
182 //- / ( msol->GetOuterRadiusMinusZ() - msol->GetInnerRadiusMinusZ());
183 G4double fwidthPlus = CalculateWidth( msol->GetOuterRadiusPlusZ()
184 - msol->GetInnerRadiusPlusZ(), fnDiv, foffset );
185 G4double pRMin2 = msol->GetInnerRadiusPlusZ()
186 + foffset + fwidthPlus * copyNo;
187 G4double pRMax2 = msol->GetInnerRadiusPlusZ()
188 + foffset + fwidthPlus * (copyNo+1);
189 G4double pDz = msol->GetZHalfLength();
190
191 G4double d_half_gap = fhgap * pRMax2 / pRMax1;
192 //- already rotated double pSR = foffset + copyNo*fwidth;
193 G4double pSPhi = msol->GetStartPhiAngle();
194 G4double pDPhi = msol->GetDeltaPhiAngle();;
195
196 cons.SetInnerRadiusMinusZ( pRMin1 + fhgap );
197 cons.SetOuterRadiusMinusZ( pRMax1 - fhgap );
198 cons.SetInnerRadiusPlusZ( pRMin2 + d_half_gap );
199 cons.SetOuterRadiusPlusZ( pRMax2 - d_half_gap );
200 cons.SetZHalfLength( pDz );
201 cons.SetStartPhiAngle( pSPhi, false );
202 cons.SetDeltaPhiAngle( pDPhi );
203
204#ifdef G4DIVDEBUG
205 if( verbose >= 2 )
206 {
207 G4cout << " G4ParameterisationConsRho::ComputeDimensions()" << G4endl
208 << " pRMin: " << pRMin1 << " - pRMax: " << pRMax1 << G4endl;
209 if( verbose >= 4 ) cons.DumpInfo();
210 }
211#endif
212}
213
214//--------------------------------------------------------------------------
217 G4double width, G4double offset,
218 G4VSolid* msolid, DivisionType divType )
219 : G4VParameterisationCons( axis, nDiv, width, offset, msolid, divType )
220{
222 SetType( "DivisionConsPhi" );
223
224 G4Cons* msol = (G4Cons*)(fmotherSolid);
225 if( divType == DivWIDTH )
226 {
227 fnDiv = CalculateNDiv( msol->GetDeltaPhiAngle(), width, offset );
228 }
229 else if( divType == DivNDIV )
230 {
231 fwidth = CalculateWidth( msol->GetDeltaPhiAngle(), nDiv, offset );
232 }
233
234#ifdef G4DIVDEBUG
235 if( verbose >= 1 )
236 {
237 G4cout << " G4ParameterisationConsPhi no divisions " << fnDiv << " = "
238 << nDiv << G4endl
239 << " Offset " << foffset << " = " << offset << G4endl
240 << " Width " << fwidth << " = " << width << G4endl;
241 }
242#endif
243}
244
245//--------------------------------------------------------------------------
247{
248}
249
250//------------------------------------------------------------------------
252{
253 G4Cons* msol = (G4Cons*)(fmotherSolid);
254 return msol->GetDeltaPhiAngle();
255}
256
257//--------------------------------------------------------------------------
258void
260ComputeTransformation( const G4int copyNo, G4VPhysicalVolume *physVol ) const
261{
262 //----- translation
263 G4ThreeVector origin(0.,0.,0.);
264 //----- set translation
265 physVol->SetTranslation( origin );
266
267 //----- calculate rotation matrix (so that all volumes point to the centre)
268 G4double posi = foffset + copyNo*fwidth;
269
270#ifdef G4DIVDEBUG
271 if( verbose >= 2 )
272 {
273 G4cout << " G4ParameterisationConsPhi - position: " << posi/deg << G4endl
274 << " Origin: " << origin << " copyNo: " << copyNo
275 << " - foffset: " << foffset/deg
276 << " - fwidth: " << fwidth/deg << G4endl
277 << " - Axis: " << faxis << G4endl;
278 }
279#endif
280
281 ChangeRotMatrix( physVol, -posi );
282}
283
284//--------------------------------------------------------------------------
285void
287ComputeDimensions( G4Cons& cons, const G4int,
288 const G4VPhysicalVolume* ) const
289{
290 G4Cons* msol = (G4Cons*)(fmotherSolid);
291
292 G4double pRMin1 = msol->GetInnerRadiusMinusZ();
293 G4double pRMax1 = msol->GetOuterRadiusMinusZ();
294 G4double pRMin2 = msol->GetInnerRadiusPlusZ();
295 G4double pRMax2 = msol->GetOuterRadiusPlusZ();
296 G4double pDz = msol->GetZHalfLength();
297
298 //- already rotated double pSPhi = foffset + copyNo*fwidth;
299 G4double pSPhi = foffset + msol->GetStartPhiAngle() + fhgap;
300 G4double pDPhi = fwidth - 2.*fhgap;
301
302 cons.SetInnerRadiusMinusZ( pRMin1 );
303 cons.SetOuterRadiusMinusZ( pRMax1 );
304 cons.SetInnerRadiusPlusZ( pRMin2 );
305 cons.SetOuterRadiusPlusZ( pRMax2 );
306 cons.SetZHalfLength( pDz );
307 cons.SetStartPhiAngle( pSPhi, false );
308 cons.SetDeltaPhiAngle( pDPhi );
309
310#ifdef G4DIVDEBUG
311 if( verbose >= 2 )
312 {
313 G4cout << " G4ParameterisationConsPhi::ComputeDimensions" << G4endl
314 << " pSPhi: " << pSPhi << " - pDPhi: " << pDPhi << G4endl;
315 if( verbose >= 4 ) cons.DumpInfo();
316 }
317#endif
318}
319
320//--------------------------------------------------------------------------
323 G4double width, G4double offset,
324 G4VSolid* msolid, DivisionType divType )
325 : G4VParameterisationCons( axis, nDiv, width, offset, msolid, divType )
326{
328 SetType( "DivisionConsZ" );
329
330 G4Cons* msol = (G4Cons*)(fmotherSolid);
331 if( divType == DivWIDTH )
332 {
333 fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(), width, offset );
334 }
335 else if( divType == DivNDIV )
336 {
337 fwidth = CalculateWidth( 2*msol->GetZHalfLength(), nDiv, offset );
338 }
339
340#ifdef G4DIVDEBUG
341 if( verbose >= 1 )
342 {
343 G4cout << " G4ParameterisationConsZ: # divisions " << fnDiv << " = "
344 << nDiv << G4endl
345 << " Offset " << foffset << " = " << offset << G4endl
346 << " Width " << fwidth << " = " << width << G4endl
347 << " - Axis: " << faxis << G4endl;
348 }
349#endif
350}
351
352//--------------------------------------------------------------------------
354{
355}
356
357//------------------------------------------------------------------------
359{
360 G4Cons* msol = (G4Cons*)(fmotherSolid);
361 return 2*msol->GetZHalfLength();
362}
363
364//--------------------------------------------------------------------------
365void
367ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol ) const
368{
369 //----- set translation: along Z axis
370 G4Cons* motherCons = (G4Cons*)(GetMotherSolid());
371 G4double posi = - motherCons->GetZHalfLength() + OffsetZ()
372 + fwidth/2 + copyNo*fwidth;
373 G4ThreeVector origin(0.,0.,posi);
374 physVol->SetTranslation( origin );
375
376 //----- calculate rotation matrix: unit
377
378#ifdef G4DIVDEBUG
379 if( verbose >= 2 )
380 {
381 G4cout << " G4ParameterisationConsZ::ComputeTransformation()" << G4endl
382 << " Origin: " << origin << " - copyNo: " << copyNo << G4endl
383 << " foffset: " << foffset << " - fwidth: " << fwidth
384 << G4endl;
385 }
386#endif
387
388 ChangeRotMatrix( physVol );
389}
390
391
392//--------------------------------------------------------------------------
393void
395ComputeDimensions( G4Cons& cons, const G4int copyNo,
396 const G4VPhysicalVolume* ) const
397{
398 G4Cons* msol = (G4Cons*)(fmotherSolid);
399
400 G4double mHalfLength = msol->GetZHalfLength() - fhgap;
401 G4double aRInner = (msol->GetInnerRadiusPlusZ()
402 - msol->GetInnerRadiusMinusZ()) / (2*mHalfLength);
403 G4double bRInner = (msol->GetInnerRadiusPlusZ()
404 + msol->GetInnerRadiusMinusZ()) / 2;
405 G4double aROuter = (msol->GetOuterRadiusPlusZ()
406 - msol->GetOuterRadiusMinusZ()) / (2*mHalfLength);
407 G4double bROuter = (msol->GetOuterRadiusPlusZ()
408 + msol->GetOuterRadiusMinusZ()) / 2;
409 G4double xMinusZ = -mHalfLength + OffsetZ() + fwidth*copyNo + fhgap;
410 G4double xPlusZ = -mHalfLength + OffsetZ() + fwidth*(copyNo+1) - fhgap;
411 cons.SetInnerRadiusMinusZ( aRInner * xMinusZ + bRInner );
412 cons.SetOuterRadiusMinusZ( aROuter * xMinusZ + bROuter );
413 cons.SetInnerRadiusPlusZ( aRInner * xPlusZ + bRInner );
414 cons.SetOuterRadiusPlusZ( aROuter * xPlusZ + bROuter );
415
416 G4double pDz = fwidth / 2. - fhgap;
417 G4double pSPhi = msol->GetStartPhiAngle();
418 G4double pDPhi = msol->GetDeltaPhiAngle();
419
420 cons.SetZHalfLength( pDz );
421 cons.SetStartPhiAngle( pSPhi, false );
422 cons.SetDeltaPhiAngle( pDPhi );
423
424#ifdef G4DIVDEBUG
425 if( verbose >= 2 )
426 {
427 G4cout << " G4ParameterisationConsZ::ComputeDimensions()" << G4endl
428 << " pDz: " << pDz << G4endl;
429 if( verbose >= 4 ) cons.DumpInfo();
430 }
431#endif
432
433}
@ JustWarning
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
Definition: G4Cons.hh:75
G4double GetOuterRadiusPlusZ() const
void SetInnerRadiusPlusZ(G4double Rmin2)
void SetZHalfLength(G4double newDz)
G4double GetStartPhiAngle() const
G4double GetDeltaPhiAngle() const
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
void SetOuterRadiusMinusZ(G4double Rmax1)
G4double GetInnerRadiusMinusZ() const
void SetOuterRadiusPlusZ(G4double Rmax2)
G4double GetInnerRadiusPlusZ() const
void SetDeltaPhiAngle(G4double newDPhi)
void SetInnerRadiusMinusZ(G4double Rmin1)
G4double GetOuterRadiusMinusZ() const
G4double GetZHalfLength() const
void ComputeDimensions(G4Cons &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const
G4ParameterisationConsPhi(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
G4ParameterisationConsRho(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
void ComputeDimensions(G4Cons &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
G4ParameterisationConsZ(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
void ComputeDimensions(G4Cons &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.) const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4VSolid * GetMotherSolid() const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
G4VParameterisationCons(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
void SetTranslation(const G4ThreeVector &v)
G4String GetName() const
void DumpInfo() const
virtual G4GeometryType GetEntityType() const =0
EAxis
Definition: geomdefs.hh:54
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41