ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ParameterisationPara.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file 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 // class G4ParameterisationPara[X,Y,Z] implementation
27 //
28 // 26.05.03 - P.Arce, Initial version
29 // 08.04.04 - I.Hrivnacova, Implemented reflection
30 // 21.04.10 - M.Asai, Added gaps
31 // --------------------------------------------------------------------
32 
34 
35 #include <iomanip>
36 
37 #include "G4PhysicalConstants.hh"
38 #include "G4ThreeVector.hh"
39 #include "G4Transform3D.hh"
40 #include "G4RotationMatrix.hh"
41 #include "G4VPhysicalVolume.hh"
42 #include "G4ReflectedSolid.hh"
43 #include "G4Para.hh"
44 
45 //--------------------------------------------------------------------------
48  G4double offset, G4VSolid* msolid,
49  DivisionType divType )
50  : G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid )
51 {
52  G4Para* msol = (G4Para*)(msolid);
53  if (msolid->GetEntityType() == "G4ReflectedSolid")
54  {
55  // Get constituent solid
56  G4VSolid* mConstituentSolid
57  = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
58  msol = (G4Para*)(mConstituentSolid);
59  fmotherSolid = msol;
60 
61  // Create a new solid with inversed parameters
62  G4Para* newSolid
63  = new G4Para(msol->GetName(),
64  msol->GetXHalfLength(),
65  msol->GetYHalfLength(),
66  msol->GetZHalfLength(),
67  std::atan(msol->GetTanAlpha()),
68  pi - msol->GetSymAxis().theta(),
69  msol->GetSymAxis().phi());
70 
71  msol = newSolid;
72  fmotherSolid = newSolid;
73  fReflectedSolid = true;
74  fDeleteSolid = true;
75  }
76 }
77 
78 //------------------------------------------------------------------------
80 {
81 }
82 
83 //------------------------------------------------------------------------
86  G4double width, G4double offset,
87  G4VSolid* msolid, DivisionType divType )
88  : G4VParameterisationPara( axis, nDiv, width, offset, msolid, divType )
89 {
91  SetType( "DivisionParaX" );
92 
93  G4Para* mpara = (G4Para*)(fmotherSolid);
94  if( divType == DivWIDTH )
95  {
96  fnDiv = CalculateNDiv( 2*mpara->GetXHalfLength(), width, offset );
97  }
98  else if( divType == DivNDIV )
99  {
100  fwidth = CalculateWidth( 2*mpara->GetXHalfLength(), nDiv, offset );
101  }
102 
103 #ifdef G4DIVDEBUG
104  if( verbose >= 1 )
105  {
106  G4cout << " G4ParameterisationParaX - # divisions " << fnDiv
107  << " = " << nDiv << G4endl
108  << " Offset " << foffset << " = " << offset << G4endl
109  << " Width " << fwidth << " = " << width << G4endl;
110  }
111 #endif
112 }
113 
114 //------------------------------------------------------------------------
116 {
117  G4Para* msol = (G4Para*)(fmotherSolid);
118  return 2*msol->GetXHalfLength();
119 }
120 
121 //------------------------------------------------------------------------
123 {
124 }
125 
126 //------------------------------------------------------------------------
127 void
129 ComputeTransformation( const G4int copyNo, G4VPhysicalVolume *physVol ) const
130 {
131  G4Para* msol = (G4Para*)(fmotherSolid );
132  G4double mdx = msol->GetXHalfLength( );
133 
134  //----- translation
135  G4ThreeVector origin(0.,0.,0.);
136  G4double posi = -mdx + foffset+(copyNo+0.5)*fwidth;
137  origin.setX( posi );
138 
139 #ifdef G4DIVDEBUG
140  if( verbose >= 2 )
141  {
142  G4cout << std::setprecision(8) << " G4ParameterisationParaX "
143  << copyNo << G4endl
144  << " Position: " << origin << " - Axis: " << faxis << G4endl;
145  }
146 #endif
147 
148  //----- set translation
149  physVol->SetTranslation( origin );
150 }
151 
152 //--------------------------------------------------------------------------
153 void
156  const G4VPhysicalVolume*) const
157 {
158  //---- The division along X of a Para will result a Para
159  G4Para* msol = (G4Para*)(fmotherSolid);
160 
161  //---- Get
162  G4double pDx = fwidth/2. - fhgap;
163  G4double pDy = msol->GetYHalfLength();
164  G4double pDz = msol->GetZHalfLength();
165  G4double pAlpha = std::atan(msol->GetTanAlpha());
166  G4double pTheta = msol->GetSymAxis().theta();
167  G4double pPhi = msol->GetSymAxis().phi();
168 
169  para.SetAllParameters ( pDx, pDy, pDz, pAlpha, pTheta, pPhi );
170 
171 #ifdef G4DIVDEBUG
172  if( verbose >= 1 )
173  {
174  G4cout << " G4ParameterisationParaX::ComputeDimensions()"
175  << " - Mother PARA " << G4endl;
176  msol->DumpInfo();
177  G4cout << " - Parameterised PARA: " << G4endl;
178  para.DumpInfo();
179  }
180 #endif
181 }
182 
183 //------------------------------------------------------------------------
186  G4double width, G4double offset,
187  G4VSolid* msolid, DivisionType divType )
188  : G4VParameterisationPara( axis, nDiv, width, offset, msolid, divType )
189 {
191  SetType( "DivisionParaY" );
192 
193  G4Para* mpara = (G4Para*)(fmotherSolid);
194  if( divType == DivWIDTH )
195  {
196  fnDiv = CalculateNDiv( 2*mpara->GetYHalfLength(), width, offset );
197  }
198  else if( divType == DivNDIV )
199  {
200  fwidth = CalculateWidth( 2*mpara->GetYHalfLength(), nDiv, offset );
201  }
202 
203 #ifdef G4DIVDEBUG
204  if( verbose >= 1 )
205  {
206  G4cout << " G4ParameterisationParaY - # divisions " << fnDiv
207  << " = " << nDiv << G4endl
208  << " Offset " << foffset << " = " << offset << G4endl
209  << " Width " << fwidth << " = " << width << G4endl;
210  }
211 #endif
212 }
213 
214 //------------------------------------------------------------------------
216 {
217 }
218 
219 //------------------------------------------------------------------------
221 {
222  G4Para* msol = (G4Para*)(fmotherSolid);
223  return 2*msol->GetYHalfLength();
224 }
225 
226 //------------------------------------------------------------------------
227 void
229 ComputeTransformation( const G4int copyNo, G4VPhysicalVolume *physVol ) const
230 {
231  G4Para* msol = (G4Para*)(fmotherSolid );
232  G4double mdy = msol->GetYHalfLength( );
233 
234  //----- translation
235  G4ThreeVector origin(0.,0.,0.);
236  G4double posiY = -mdy + foffset+(copyNo+0.5)*fwidth;
237  origin.setY( posiY );
238  G4double posiX = posiY * msol->GetTanAlpha();
239  origin.setX( posiX );
240 
241 #ifdef G4DIVDEBUG
242  if( verbose >= 2 )
243  {
244  G4cout << std::setprecision(8) << " G4ParameterisationParaY "
245  << copyNo << G4endl
246  << " Position: " << origin << " - Axis: " << faxis << G4endl;
247  }
248 #endif
249 
250  //----- set translation
251  physVol->SetTranslation( origin );
252 }
253 
254 //--------------------------------------------------------------------------
255 void
258  const G4VPhysicalVolume*) const
259 {
260  //---- The division along Y of a Para will result a Para
261  G4Para* msol = (G4Para*)(fmotherSolid);
262 
263  //---- Get
264  G4double pDx = msol->GetXHalfLength();
265  G4double pDy = fwidth/2. - fhgap;
266  G4double pDz = msol->GetZHalfLength();
267  G4double pAlpha = std::atan(msol->GetTanAlpha());
268  G4double pTheta = msol->GetSymAxis().theta();
269  G4double pPhi = msol->GetSymAxis().phi();
270 
271  para.SetAllParameters ( pDx, pDy, pDz, pAlpha, pTheta, pPhi );
272 
273 #ifdef G4DIVDEBUG
274  if( verbose >= -1 )
275  {
276  G4cout << " G4ParameterisationParaY::ComputeDimensions()"
277  << " - Mother PARA " << G4endl;
278  msol->DumpInfo();
279  G4cout << " - Parameterised PARA: " << G4endl;
280  para.DumpInfo();
281  }
282 #endif
283 }
284 
285 //------------------------------------------------------------------------
288  G4double width, G4double offset,
289  G4VSolid* msolid, DivisionType divType )
290  : G4VParameterisationPara( axis, nDiv, width, offset, msolid, divType )
291 {
293  SetType( "DivisionParaZ" );
294 
295  G4Para* mpara = (G4Para*)(fmotherSolid);
296  if( divType == DivWIDTH )
297  {
298  fnDiv = CalculateNDiv( 2*mpara->GetZHalfLength(), width, offset );
299  }
300  else if( divType == DivNDIV )
301  {
302  fwidth = CalculateWidth( 2*mpara->GetZHalfLength(), nDiv, offset );
303  }
304 
305 #ifdef G4DIVDEBUG
306  if( verbose >= -1 )
307  {
308  G4cout << " G4ParameterisationParaZ - # divisions " << fnDiv
309  << " = " << nDiv << G4endl
310  << " Offset " << foffset << " = " << offset << G4endl
311  << " Width " << fwidth << " = " << width << G4endl;
312  }
313 #endif
314 }
315 
316 //------------------------------------------------------------------------
318 {
319 }
320 
321 //------------------------------------------------------------------------
323 {
324  G4Para* msol = (G4Para*)(fmotherSolid);
325  return 2*msol->GetZHalfLength();
326 }
327 
328 //------------------------------------------------------------------------
329 void
331 ComputeTransformation( const G4int copyNo, G4VPhysicalVolume *physVol ) const
332 {
333  G4Para* msol = (G4Para*)(fmotherSolid );
334  G4double mdz = msol->GetZHalfLength( );
335 
336  //----- translation
337  G4double posi = -mdz + OffsetZ() + (copyNo+0.5)*fwidth;
338  G4ThreeVector symAxis = msol->GetSymAxis();
339  G4ThreeVector origin( symAxis * posi / symAxis.z() );
340 
341 #ifdef G4DIVDEBUG
342  if( verbose >= 2 )
343  {
344  G4cout << std::setprecision(8) << " G4ParameterisationParaZ "
345  << copyNo << G4endl
346  << " Position: " << origin << " - Axis: " << faxis << G4endl;
347  }
348 #endif
349 
350  //----- set translation
351  physVol->SetTranslation( origin );
352 }
353 
354 //--------------------------------------------------------------------------
355 void
358  const G4VPhysicalVolume*) const
359 {
360  //---- The division along Z of a Para will result a Para
361  G4Para* msol = (G4Para*)(fmotherSolid);
362 
363  //---- Get
364  G4double pDx = msol->GetXHalfLength();
365  G4double pDy = msol->GetYHalfLength();
366  G4double pDz = fwidth/2. - fhgap;
367  G4double pAlpha = std::atan(msol->GetTanAlpha());
368  G4double pTheta = msol->GetSymAxis().theta();
369  G4double pPhi = msol->GetSymAxis().phi();
370 
371  para.SetAllParameters ( pDx, pDy, pDz, pAlpha, pTheta, pPhi );
372 
373 #ifdef G4DIVDEBUG
374  if( verbose >= -1 )
375  {
376  G4cout << " G4ParameterisationParaZ::ComputeDimensions()"
377  << " - Mother PARA " << G4endl;
378  msol->DumpInfo();
379  G4cout << " - Parameterised PARA: " << G4endl;
380  para.DumpInfo();
381  }
382 #endif
383 }