Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4ParameterisationPara.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: G4ParameterisationPara.cc 66356 2012-12-18 09:02:32Z gcosmo $
28 //
29 // class G4ParameterisationPara 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 
40 #include "G4PhysicalConstants.hh"
41 #include "G4ThreeVector.hh"
42 #include "G4Transform3D.hh"
43 #include "G4RotationMatrix.hh"
44 #include "G4VPhysicalVolume.hh"
45 #include "G4ReflectedSolid.hh"
46 #include "G4Para.hh"
47 
48 //--------------------------------------------------------------------------
51  G4double offset, G4VSolid* msolid,
52  DivisionType divType )
53  : G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid )
54 {
55  G4Para* msol = (G4Para*)(msolid);
56  if (msolid->GetEntityType() == "G4ReflectedSolid")
57  {
58  // Get constituent solid
59  G4VSolid* mConstituentSolid
60  = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
61  msol = (G4Para*)(mConstituentSolid);
62  fmotherSolid = msol;
63 
64  // Create a new solid with inversed parameters
65  G4Para* newSolid
66  = new G4Para(msol->GetName(),
67  msol->GetXHalfLength(),
68  msol->GetYHalfLength(),
69  msol->GetZHalfLength(),
70  std::atan(msol->GetTanAlpha()),
71  pi - msol->GetSymAxis().theta(),
72  msol->GetSymAxis().phi());
73 
74  msol = newSolid;
75  fmotherSolid = newSolid;
76  fReflectedSolid = true;
77  fDeleteSolid = true;
78  }
79 }
80 
81 //------------------------------------------------------------------------
83 {
84 }
85 
86 //------------------------------------------------------------------------
89  G4double width, G4double offset,
90  G4VSolid* msolid, DivisionType divType )
91  : G4VParameterisationPara( axis, nDiv, width, offset, msolid, divType )
92 {
94  SetType( "DivisionParaX" );
95 
96  G4Para* mpara = (G4Para*)(fmotherSolid);
97  if( divType == DivWIDTH )
98  {
99  fnDiv = CalculateNDiv( 2*mpara->GetXHalfLength(), width, offset );
100  }
101  else if( divType == DivNDIV )
102  {
103  fwidth = CalculateWidth( 2*mpara->GetXHalfLength(), nDiv, offset );
104  }
105 
106 #ifdef G4DIVDEBUG
107  if( verbose >= 1 )
108  {
109  G4cout << " G4ParameterisationParaX - # divisions " << fnDiv
110  << " = " << nDiv << G4endl
111  << " Offset " << foffset << " = " << offset << G4endl
112  << " Width " << fwidth << " = " << width << G4endl;
113  }
114 #endif
115 }
116 
117 //------------------------------------------------------------------------
119 {
120  G4Para* msol = (G4Para*)(fmotherSolid);
121  return 2*msol->GetXHalfLength();
122 }
123 
124 //------------------------------------------------------------------------
126 {
127 }
128 
129 //------------------------------------------------------------------------
130 void
132 ComputeTransformation( const G4int copyNo, G4VPhysicalVolume *physVol ) const
133 {
134  G4Para* msol = (G4Para*)(fmotherSolid );
135  G4double mdx = msol->GetXHalfLength( );
136 
137  //----- translation
138  G4ThreeVector origin(0.,0.,0.);
139  G4double posi = -mdx + foffset+(copyNo+0.5)*fwidth;
140  origin.setX( posi );
141 
142 #ifdef G4DIVDEBUG
143  if( verbose >= 2 )
144  {
145  G4cout << std::setprecision(8) << " G4ParameterisationParaX "
146  << copyNo << G4endl
147  << " Position: " << origin << " - Axis: " << faxis << G4endl;
148  }
149 #endif
150 
151  //----- set translation
152  physVol->SetTranslation( origin );
153 }
154 
155 //--------------------------------------------------------------------------
156 void
159  const G4VPhysicalVolume*) const
160 {
161  //---- The division along X of a Para will result a Para
162  G4Para* msol = (G4Para*)(fmotherSolid);
163 
164  //---- Get
165  G4double pDx = fwidth/2. - fhgap;
166  G4double pDy = msol->GetYHalfLength();
167  G4double pDz = msol->GetZHalfLength();
168  G4double pAlpha = std::atan(msol->GetTanAlpha());
169  G4double pTheta = msol->GetSymAxis().theta();
170  G4double pPhi = msol->GetSymAxis().phi();
171 
172  para.SetAllParameters ( pDx, pDy, pDz, pAlpha, pTheta, pPhi );
173 
174 #ifdef G4DIVDEBUG
175  if( verbose >= 1 )
176  {
177  G4cout << " G4ParameterisationParaX::ComputeDimensions()"
178  << " - Mother PARA " << G4endl;
179  msol->DumpInfo();
180  G4cout << " - Parameterised PARA: " << G4endl;
181  para.DumpInfo();
182  }
183 #endif
184 }
185 
186 //------------------------------------------------------------------------
189  G4double width, G4double offset,
190  G4VSolid* msolid, DivisionType divType )
191  : G4VParameterisationPara( axis, nDiv, width, offset, msolid, divType )
192 {
194  SetType( "DivisionParaY" );
195 
196  G4Para* mpara = (G4Para*)(fmotherSolid);
197  if( divType == DivWIDTH )
198  {
199  fnDiv = CalculateNDiv( 2*mpara->GetYHalfLength(), width, offset );
200  }
201  else if( divType == DivNDIV )
202  {
203  fwidth = CalculateWidth( 2*mpara->GetYHalfLength(), nDiv, offset );
204  }
205 
206 #ifdef G4DIVDEBUG
207  if( verbose >= 1 )
208  {
209  G4cout << " G4ParameterisationParaY - # divisions " << fnDiv
210  << " = " << nDiv << G4endl
211  << " Offset " << foffset << " = " << offset << G4endl
212  << " Width " << fwidth << " = " << width << G4endl;
213  }
214 #endif
215 }
216 
217 //------------------------------------------------------------------------
219 {
220 }
221 
222 //------------------------------------------------------------------------
224 {
225  G4Para* msol = (G4Para*)(fmotherSolid);
226  return 2*msol->GetYHalfLength();
227 }
228 
229 //------------------------------------------------------------------------
230 void
232 ComputeTransformation( const G4int copyNo, G4VPhysicalVolume *physVol ) const
233 {
234  G4Para* msol = (G4Para*)(fmotherSolid );
235  G4double mdy = msol->GetYHalfLength( );
236 
237  //----- translation
238  G4ThreeVector origin(0.,0.,0.);
239  G4double posiY = -mdy + foffset+(copyNo+0.5)*fwidth;
240  origin.setY( posiY );
241  G4double posiX = posiY * msol->GetTanAlpha();
242  origin.setX( posiX );
243 
244 #ifdef G4DIVDEBUG
245  if( verbose >= 2 )
246  {
247  G4cout << std::setprecision(8) << " G4ParameterisationParaY "
248  << copyNo << G4endl
249  << " Position: " << origin << " - Axis: " << faxis << G4endl;
250  }
251 #endif
252 
253  //----- set translation
254  physVol->SetTranslation( origin );
255 }
256 
257 //--------------------------------------------------------------------------
258 void
261  const G4VPhysicalVolume*) const
262 {
263  //---- The division along Y of a Para will result a Para
264  G4Para* msol = (G4Para*)(fmotherSolid);
265 
266  //---- Get
267  G4double pDx = msol->GetXHalfLength();
268  G4double pDy = fwidth/2. - fhgap;
269  G4double pDz = msol->GetZHalfLength();
270  G4double pAlpha = std::atan(msol->GetTanAlpha());
271  G4double pTheta = msol->GetSymAxis().theta();
272  G4double pPhi = msol->GetSymAxis().phi();
273 
274  para.SetAllParameters ( pDx, pDy, pDz, pAlpha, pTheta, pPhi );
275 
276 #ifdef G4DIVDEBUG
277  if( verbose >= -1 )
278  {
279  G4cout << " G4ParameterisationParaY::ComputeDimensions()"
280  << " - Mother PARA " << G4endl;
281  msol->DumpInfo();
282  G4cout << " - Parameterised PARA: " << G4endl;
283  para.DumpInfo();
284  }
285 #endif
286 }
287 
288 //------------------------------------------------------------------------
291  G4double width, G4double offset,
292  G4VSolid* msolid, DivisionType divType )
293  : G4VParameterisationPara( axis, nDiv, width, offset, msolid, divType )
294 {
296  SetType( "DivisionParaZ" );
297 
298  G4Para* mpara = (G4Para*)(fmotherSolid);
299  if( divType == DivWIDTH )
300  {
301  fnDiv = CalculateNDiv( 2*mpara->GetZHalfLength(), width, offset );
302  }
303  else if( divType == DivNDIV )
304  {
305  fwidth = CalculateWidth( 2*mpara->GetZHalfLength(), nDiv, offset );
306  }
307 
308 #ifdef G4DIVDEBUG
309  if( verbose >= -1 )
310  {
311  G4cout << " G4ParameterisationParaZ - # divisions " << fnDiv
312  << " = " << nDiv << G4endl
313  << " Offset " << foffset << " = " << offset << G4endl
314  << " Width " << fwidth << " = " << width << G4endl;
315  }
316 #endif
317 }
318 
319 //------------------------------------------------------------------------
321 {
322 }
323 
324 //------------------------------------------------------------------------
326 {
327  G4Para* msol = (G4Para*)(fmotherSolid);
328  return 2*msol->GetZHalfLength();
329 }
330 
331 //------------------------------------------------------------------------
332 void
334 ComputeTransformation( const G4int copyNo, G4VPhysicalVolume *physVol ) const
335 {
336  G4Para* msol = (G4Para*)(fmotherSolid );
337  G4double mdz = msol->GetZHalfLength( );
338 
339  //----- translation
340  G4double posi = -mdz + OffsetZ() + (copyNo+0.5)*fwidth;
341  G4ThreeVector symAxis = msol->GetSymAxis();
342  G4ThreeVector origin( symAxis * posi / symAxis.z() );
343 
344 #ifdef G4DIVDEBUG
345  if( verbose >= 2 )
346  {
347  G4cout << std::setprecision(8) << " G4ParameterisationParaZ "
348  << copyNo << G4endl
349  << " Position: " << origin << " - Axis: " << faxis << G4endl;
350  }
351 #endif
352 
353  //----- set translation
354  physVol->SetTranslation( origin );
355 }
356 
357 //--------------------------------------------------------------------------
358 void
361  const G4VPhysicalVolume*) const
362 {
363  //---- The division along Z of a Para will result a Para
364  G4Para* msol = (G4Para*)(fmotherSolid);
365 
366  //---- Get
367  G4double pDx = msol->GetXHalfLength();
368  G4double pDy = msol->GetYHalfLength();
369  G4double pDz = fwidth/2. - fhgap;
370  G4double pAlpha = std::atan(msol->GetTanAlpha());
371  G4double pTheta = msol->GetSymAxis().theta();
372  G4double pPhi = msol->GetSymAxis().phi();
373 
374  para.SetAllParameters ( pDx, pDy, pDz, pAlpha, pTheta, pPhi );
375 
376 #ifdef G4DIVDEBUG
377  if( verbose >= -1 )
378  {
379  G4cout << " G4ParameterisationParaZ::ComputeDimensions()"
380  << " - Mother PARA " << G4endl;
381  msol->DumpInfo();
382  G4cout << " - Parameterised PARA: " << G4endl;
383  para.DumpInfo();
384  }
385 #endif
386 }
G4ParameterisationParaY(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
G4double GetTanAlpha() const
void ComputeDimensions(G4Para &para, const G4int copyNo, const G4VPhysicalVolume *pv) const
#define G4endl
Definition: G4ios.hh:61
G4double GetXHalfLength() const
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
double z() const
void SetTranslation(const G4ThreeVector &v)
G4ThreeVector GetSymAxis() const
void SetAllParameters(G4double pDx, G4double pDy, G4double pDz, G4double pAlpha, G4double pTheta, G4double pPhi)
Definition: G4Para.cc:200
void setX(double)
double theta() const
double G4double
Definition: G4Types.hh:76
#define width
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
G4VParameterisationPara(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
void ComputeDimensions(G4Para &para, const G4int copyNo, const G4VPhysicalVolume *pv) const
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
int G4int
Definition: G4Types.hh:78
double phi() const
EAxis
Definition: geomdefs.hh:54
void ComputeDimensions(G4Para &para, const G4int copyNo, const G4VPhysicalVolume *pv) const
G4String GetName() const
void DumpInfo() const
G4GLOB_DLL std::ostream G4cout
G4ParameterisationParaX(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4double GetZHalfLength() const
static constexpr double pi
Definition: G4SIunits.hh:75
G4ParameterisationParaZ(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
virtual G4GeometryType GetEntityType() const =0
void SetType(const G4String &type)
Definition: G4Para.hh:86
void setY(double)
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
G4double GetYHalfLength() const