Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4ParameterisationTubs.cc
이 파일의 문서화 페이지로 가기
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: G4ParameterisationTubs.cc 66356 2012-12-18 09:02:32Z gcosmo $
28 //
29 // class G4ParameterisationTubs 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 "G4Tubs.hh"
45 
46 //--------------------------------------------------------------------------
49  G4double offset, G4VSolid* msolid,
50  DivisionType divType )
51  : G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid )
52 {
53  G4Tubs* msol = (G4Tubs*)(msolid);
54  if (msolid->GetEntityType() == "G4ReflectedSolid")
55  {
56  //----- get constituent solid
57  G4VSolid* mConstituentSolid
58  = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
59  msol = (G4Tubs*)(mConstituentSolid);
60  fmotherSolid = msol;
61  fReflectedSolid = true;
62  }
63 }
64 
65 //------------------------------------------------------------------------
67 {
68 }
69 
70 //--------------------------------------------------------------------------
73  G4double width, G4double offset,
74  G4VSolid* msolid, DivisionType divType )
75  : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
76 {
78  SetType( "DivisionTubsRho" );
79 
80  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
81  if( divType == DivWIDTH )
82  {
83  fnDiv = CalculateNDiv( msol->GetOuterRadius() - msol->GetInnerRadius(),
84  width, offset );
85  }
86  else if( divType == DivNDIV )
87  {
89  nDiv, offset );
90  }
91 
92 #ifdef G4DIVDEBUG
93  if( verbose >= 1 )
94  {
95  G4cout << " G4ParameterisationTubsRho - no divisions " << fnDiv << " = "
96  << nDiv << G4endl
97  << " Offset " << foffset << " = " << offset << G4endl
98  << " Width " << fwidth << " = " << width << G4endl
99  << " DivType " << divType << G4endl;
100  }
101 #endif
102 }
103 
104 //--------------------------------------------------------------------------
106 {
107 }
108 
109 //------------------------------------------------------------------------
111 {
112  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
113  return msol->GetOuterRadius() - msol->GetInnerRadius();
114 }
115 
116 
117 //--------------------------------------------------------------------------
118 void
121 {
122  //----- translation
123  G4ThreeVector origin(0.,0.,0.);
124  //----- set translation
125  physVol->SetTranslation( origin );
126 
127  //----- calculate rotation matrix: unit
128 
129 #ifdef G4DIVDEBUG
130  if( verbose >= 2 )
131  {
132  G4cout << " G4ParameterisationTubsRho " << G4endl
133  << " Offset: " << foffset/deg
134  << " - Width: " << fwidth/deg << G4endl;
135  }
136 #endif
137 
138  ChangeRotMatrix( physVol );
139 
140 #ifdef G4DIVDEBUG
141  if( verbose >= 2 )
142  {
143  G4cout << std::setprecision(8) << " G4ParameterisationTubsRho " << G4endl
144  << " Position: " << origin << " - Width: " << fwidth
145  << " - Axis " << faxis << G4endl;
146  }
147 #endif
148 }
149 
150 //--------------------------------------------------------------------------
151 void
153 ComputeDimensions( G4Tubs& tubs, const G4int copyNo,
154  const G4VPhysicalVolume* ) const
155 {
156  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
157 
158  G4double pRMin = msol->GetInnerRadius() + foffset + fwidth * copyNo + fhgap;
159  G4double pRMax = msol->GetInnerRadius() + foffset + fwidth * (copyNo+1) - fhgap;
160  G4double pDz = msol->GetZHalfLength();
161  //- already rotated G4double pSR = foffset + copyNo*fwidth;
162  G4double pSPhi = msol->GetStartPhiAngle();
163  G4double pDPhi = msol->GetDeltaPhiAngle();;
164 
165  tubs.SetInnerRadius( pRMin );
166  tubs.SetOuterRadius( pRMax );
167  tubs.SetZHalfLength( pDz );
168  tubs.SetStartPhiAngle( pSPhi, false );
169  tubs.SetDeltaPhiAngle( pDPhi );
170 
171 #ifdef G4DIVDEBUG
172  if( verbose >= 2 )
173  {
174  G4cout << " G4ParameterisationTubsRho::ComputeDimensions()" << G4endl
175  << " pRMin: " << pRMin << " - pRMax: " << pRMax << G4endl;
176  tubs.DumpInfo();
177  }
178 #endif
179 }
180 
181 //--------------------------------------------------------------------------
184  G4double width, G4double offset,
185  G4VSolid* msolid, DivisionType divType )
186  : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
187 {
189  SetType( "DivisionTubsPhi" );
190 
191  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
192  if( divType == DivWIDTH )
193  {
194  fnDiv = CalculateNDiv( msol->GetDeltaPhiAngle(), width, offset );
195  }
196  else if( divType == DivNDIV )
197  {
198  fwidth = CalculateWidth( msol->GetDeltaPhiAngle(), nDiv, offset );
199  }
200 
201 #ifdef G4DIVDEBUG
202  if( verbose >= 1 )
203  {
204  G4cout << " G4ParameterisationTubsPhi no divisions " << fnDiv << " = "
205  << nDiv << G4endl
206  << " Offset " << foffset << " = " << offset << G4endl
207  << " Width " << fwidth << " = " << width << G4endl;
208  }
209 #endif
210 }
211 
212 //--------------------------------------------------------------------------
214 {
215 }
216 
217 //------------------------------------------------------------------------
219 {
220  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
221  return msol->GetDeltaPhiAngle();
222 }
223 
224 //--------------------------------------------------------------------------
225 void
227 ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
228 {
229  //----- translation
230  G4ThreeVector origin(0.,0.,0.);
231  //----- set translation
232  physVol->SetTranslation( origin );
233 
234  //----- calculate rotation matrix (so that all volumes point to the centre)
235  G4double posi = foffset + copyNo*fwidth;
236 
237 #ifdef G4DIVDEBUG
238  if( verbose >= 2 )
239  {
240  G4cout << " G4ParameterisationTubsPhi - position: " << posi/deg << G4endl
241  << " copyNo: " << copyNo << " - foffset: " << foffset/deg
242  << " - fwidth: " << fwidth/deg << G4endl;
243  }
244 #endif
245 
246  ChangeRotMatrix( physVol, -posi );
247 
248 #ifdef G4DIVDEBUG
249  if( verbose >= 2 )
250  {
251  G4cout << std::setprecision(8) << " G4ParameterisationTubsPhi " << copyNo
252  << G4endl
253  << " Position: " << origin << " - Width: " << fwidth
254  << " - Axis: " << faxis << G4endl;
255  }
256 #endif
257 }
258 
259 //--------------------------------------------------------------------------
260 void
263  const G4VPhysicalVolume* ) const
264 {
265  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
266 
267  G4double pRMin = msol->GetInnerRadius();
268  G4double pRMax = msol->GetOuterRadius();
269  G4double pDz = msol->GetZHalfLength();
270  //----- already rotated in 'ComputeTransformation'
271  G4double pSPhi = msol->GetStartPhiAngle() + fhgap;
272  G4double pDPhi = fwidth - 2.*fhgap;
273 
274  tubs.SetInnerRadius( pRMin );
275  tubs.SetOuterRadius( pRMax );
276  tubs.SetZHalfLength( pDz );
277  tubs.SetStartPhiAngle( pSPhi, false );
278  tubs.SetDeltaPhiAngle( pDPhi );
279 
280 #ifdef G4DIVDEBUG
281  if( verbose >= 2 )
282  {
283  G4cout << " G4ParameterisationTubsPhi::ComputeDimensions" << G4endl
284  << " pSPhi: " << pSPhi << " - pDPhi: " << pDPhi << G4endl;
285  tubs.DumpInfo();
286  }
287 #endif
288 }
289 
290 //--------------------------------------------------------------------------
293  G4double width, G4double offset,
294  G4VSolid* msolid, DivisionType divType )
295  : G4VParameterisationTubs( axis, nDiv, width, offset, msolid, divType )
296 {
298  SetType( "DivisionTubsZ" );
299 
300  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
301  if( divType == DivWIDTH )
302  {
303  fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(), width, offset );
304  }
305  else if( divType == DivNDIV )
306  {
307  fwidth = CalculateWidth( 2*msol->GetZHalfLength(), nDiv, offset );
308  }
309 
310 #ifdef G4DIVDEBUG
311  if( verbose >= 1 )
312  {
313  G4cout << " G4ParameterisationTubsZ: # divisions " << fnDiv << " = "
314  << nDiv << G4endl
315  << " Offset " << foffset << " = " << offset << G4endl
316  << " Width " << fwidth << " = " << width << G4endl;
317  }
318 #endif
319 }
320 
321 //--------------------------------------------------------------------------
323 {
324 }
325 
326 //------------------------------------------------------------------------
328 {
329  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
330  return 2*msol->GetZHalfLength();
331 }
332 
333 //--------------------------------------------------------------------------
334 void
336 ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
337 {
338  //----- set translation: along Z axis
339  G4Tubs* motherTubs = (G4Tubs*)(fmotherSolid);
340  G4double posi = - motherTubs->GetZHalfLength() + OffsetZ()
341  + fwidth/2 + copyNo*fwidth;
342  G4ThreeVector origin(0.,0.,posi);
343  physVol->SetTranslation( origin );
344 
345  //----- calculate rotation matrix: unit
346 
347 #ifdef G4DIVDEBUG
348  if( verbose >= 2 )
349  {
350  G4cout << " G4ParameterisationTubsZ::ComputeTransformation()" << G4endl
351  << " Position: " << posi << " - copyNo: " << copyNo << G4endl
352  << " foffset " << foffset/deg << " - fwidth " << fwidth/deg
353  << G4endl;
354  }
355 #endif
356 
357  ChangeRotMatrix( physVol );
358 
359 #ifdef G4DIVDEBUG
360  if( verbose >= 2 )
361  {
362  G4cout << std::setprecision(8) << " G4ParameterisationTubsZ " << copyNo
363  << G4endl
364  << " Position: " << origin << " - Width: " << fwidth
365  << " - Axis: " << faxis << G4endl;
366  }
367 #endif
368 }
369 
370 //--------------------------------------------------------------------------
371 void
374  const G4VPhysicalVolume* ) const
375 {
376  G4Tubs* msol = (G4Tubs*)(fmotherSolid);
377 
378  G4double pRMin = msol->GetInnerRadius();
379  G4double pRMax = msol->GetOuterRadius();
380  // G4double pDz = msol->GetZHalfLength() / GetNoDiv();
381  G4double pDz = fwidth/2. - fhgap;
382  G4double pSPhi = msol->GetStartPhiAngle();
383  G4double pDPhi = msol->GetDeltaPhiAngle();
384 
385  tubs.SetInnerRadius( pRMin );
386  tubs.SetOuterRadius( pRMax );
387  tubs.SetZHalfLength( pDz );
388  tubs.SetStartPhiAngle( pSPhi, false );
389  tubs.SetDeltaPhiAngle( pDPhi );
390 
391 #ifdef G4DIVDEBUG
392  if( verbose >= 2 )
393  {
394  G4cout << " G4ParameterisationTubsZ::ComputeDimensions()" << G4endl
395  << " pDz: " << pDz << G4endl;
396  tubs.DumpInfo();
397  }
398 #endif
399 }
void SetInnerRadius(G4double newRMin)
void SetOuterRadius(G4double newRMax)
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
void ComputeDimensions(G4Tubs &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const
Definition: G4Tubs.hh:85
#define G4endl
Definition: G4ios.hh:61
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
void SetTranslation(const G4ThreeVector &v)
G4double GetDeltaPhiAngle() const
G4VParameterisationTubs(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
double G4double
Definition: G4Types.hh:76
#define width
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
G4ParameterisationTubsRho(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
static constexpr double deg
Definition: G4SIunits.hh:152
void SetZHalfLength(G4double newDz)
G4double GetInnerRadius() const
G4double GetZHalfLength() const
G4double GetStartPhiAngle() const
void ChangeRotMatrix(G4VPhysicalVolume *physVol, G4double rotZ=0.) const
void SetDeltaPhiAngle(G4double newDPhi)
void ComputeDimensions(G4Tubs &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const
int G4int
Definition: G4Types.hh:78
G4ParameterisationTubsPhi(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
EAxis
Definition: geomdefs.hh:54
void DumpInfo() const
G4double GetOuterRadius() const
G4GLOB_DLL std::ostream G4cout
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
void ComputeDimensions(G4Tubs &tubs, const G4int copyNo, const G4VPhysicalVolume *physVol) const
virtual G4GeometryType GetEntityType() const =0
void SetType(const G4String &type)
G4ParameterisationTubsZ(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType)
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const