Geant4  10.02
G4ParameterisationTrd.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: G4ParameterisationTrd.cc 68040 2013-03-13 14:19:04Z gcosmo $
28 //
29 // class G4ParameterisationTrd 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 
36 #include "G4ParameterisationTrd.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 "G4Trd.hh"
45 #include "G4Trap.hh"
46 
47 //--------------------------------------------------------------------------
50  G4double offset, G4VSolid* msolid,
51  DivisionType divType )
52  : G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid ),
53  bDivInTrap(false)
54 {
55  G4Trd* msol = (G4Trd*)(msolid);
56  if (msolid->GetEntityType() == "G4ReflectedSolid")
57  {
58  // Get constituent solid
59  G4VSolid* mConstituentSolid
60  = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
61  msol = (G4Trd*)(mConstituentSolid);
62 
63  // Create a new solid with inversed parameters
64  G4Trd* newSolid
65  = new G4Trd(msol->GetName(),
66  msol->GetXHalfLength2(), msol->GetXHalfLength1(),
67  msol->GetYHalfLength2(), msol->GetYHalfLength1(),
68  msol->GetZHalfLength());
69  msol = newSolid;
70  fmotherSolid = newSolid;
71  fReflectedSolid = true;
72  fDeleteSolid = true;
73  }
74 }
75 
76 //------------------------------------------------------------------------
78 {
79 }
80 
81 //------------------------------------------------------------------------
84  G4double width, G4double offset,
85  G4VSolid* msolid, DivisionType divType )
86  : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
87 {
89  SetType( "DivisionTrdX" );
90 
91  G4Trd* msol = (G4Trd*)(fmotherSolid);
92  if( divType == DivWIDTH )
93  {
95  width, offset );
96  }
97  else if( divType == DivNDIV )
98  {
100  nDiv, offset );
101  }
102 
103 #ifdef G4DIVDEBUG
104  if( verbose >= 1 )
105  {
106  G4cout << " G4ParameterisationTrdX - ## divisions " << fnDiv << " = "
107  << nDiv << G4endl
108  << " Offset " << foffset << " = " << offset << G4endl
109  << " Width " << fwidth << " = " << width << G4endl;
110  }
111 #endif
112 
113  G4double mpDx1 = msol->GetXHalfLength1();
114  G4double mpDx2 = msol->GetXHalfLength2();
115  if( std::fabs(mpDx1 - mpDx2) > kCarTolerance )
116  {
117  bDivInTrap = true;
118  }
119 }
120 
121 //------------------------------------------------------------------------
123 {
124 }
125 
126 //------------------------------------------------------------------------
128 {
129  G4Trd* msol = (G4Trd*)(fmotherSolid);
130  return (msol->GetXHalfLength1()+msol->GetXHalfLength2());
131 }
132 
133 //------------------------------------------------------------------------
134 void
137  G4VPhysicalVolume *physVol ) const
138 {
139  G4Trd* msol = (G4Trd*)(fmotherSolid );
140  G4double mdx = ( msol->GetXHalfLength1() + msol->GetXHalfLength2() ) / 2.;
141 
142  //----- translation
143  G4ThreeVector origin(0.,0.,0.);
144  G4double posi;
145  if( !bDivInTrap )
146  {
147  posi = -mdx + foffset + (copyNo+0.5)*fwidth;
148  }
149  else
150  {
151  G4double aveHL = (msol->GetXHalfLength1()+msol->GetXHalfLength2())/2.;
152  posi = - aveHL + foffset + (copyNo+0.5)*aveHL/fnDiv*2;
153  }
154  if( faxis == kXAxis )
155  {
156  origin.setX( posi );
157  }
158  else
159  {
160  std::ostringstream message;
161  message << "Only axes along X are allowed ! Axis: " << faxis;
162  G4Exception("G4ParameterisationTrdX::ComputeTransformation()",
163  "GeomDiv0002", FatalException, message);
164  }
165 
166 #ifdef G4DIVDEBUG
167  if( verbose >= 2 )
168  {
169  G4cout << std::setprecision(8)
170  << " G4ParameterisationTrdX::ComputeTransformation() "
171  << copyNo << G4endl
172  << " Position: " << origin << " - Axis: " << faxis << G4endl;
173  }
174 #endif
175 
176  //----- set translation
177  physVol->SetTranslation( origin );
178 }
179 
180 //--------------------------------------------------------------------------
181 void
183 ComputeDimensions( G4Trd& trd, const G4int, const G4VPhysicalVolume* ) const
184 {
185  G4Trd* msol = (G4Trd*)(fmotherSolid);
186 
187  G4double pDy1 = msol->GetYHalfLength1();
188  G4double pDy2 = msol->GetYHalfLength2();
189  G4double pDz = msol->GetZHalfLength();
190  G4double pDx = fwidth/2. - fhgap;
191 
192  trd.SetAllParameters ( pDx, pDx, pDy1, pDy2, pDz );
193 
194 #ifdef G4DIVDEBUG
195  if( verbose >= 2 )
196  {
197  G4cout << " G4ParameterisationTrdX::ComputeDimensions():"
198  << copyNo << G4endl;
199  trd.DumpInfo();
200  }
201 #endif
202 }
203 
206 {
207  if( bDivInTrap )
208  {
210  }
211  else
212  {
213  return fmotherSolid;
214  }
215 }
216 
217 
218 //--------------------------------------------------------------------------
219 void
221  const G4VPhysicalVolume* ) const
222 {
223  G4Trd* msol = (G4Trd*)(fmotherSolid);
224  G4double pDy1 = msol->GetYHalfLength1();
225  G4double pDy2 = msol->GetYHalfLength2();
226  G4double pDz = msol->GetZHalfLength();
227  G4double pDx1 = msol->GetXHalfLength1()/fnDiv;
228  G4double pDx2 = msol->GetXHalfLength2()/fnDiv;
229 
230  G4double cxy1 = -msol->GetXHalfLength1() + foffset
231  + (copyNo+0.5)*pDx1*2;// centre of the side at y=-pDy1
232  G4double cxy2 = -msol->GetXHalfLength2() + foffset
233  + (copyNo+0.5)*pDx2*2;// centre of the side at y=+pDy1
234  G4double alp = std::atan( (cxy2-cxy1)/pDz );
235 
236  trap.SetAllParameters ( pDz,
237  0.,
238  0.,
239  pDy1,
240  pDx1,
241  pDx2,
242  alp,
243  pDy2,
244  pDx1 - fhgap,
245  pDx2 - fhgap * pDx2/pDx1,
246  alp);
247 
248 #ifdef G4DIVDEBUG
249  if( verbose >= 2 )
250  {
251  G4cout << " G4ParameterisationTrdX::ComputeDimensions():"
252  << copyNo << G4endl;
253  trap.DumpInfo();
254  }
255 #endif
256 }
257 
258 //--------------------------------------------------------------------------
260 {
262 /*
263  G4Trd* msol = (G4Trd*)(fmotherSolid);
264 
265  G4double mpDx1 = msol->GetXHalfLength1();
266  G4double mpDx2 = msol->GetXHalfLength2();
267  bDivInTrap = false;
268 
269  if( std::fabs(mpDx1 - mpDx2) > kCarTolerance )
270  {
271  std::ostringstream message;
272  message << "Invalid solid specification. NOT supported." << G4endl
273  << "Making a division of a TRD along axis X," << G4endl
274  << "while the X half lengths are not equal," << G4endl
275  << "is not (yet) supported. It will result" << G4endl
276  << "in non-equal division solids.";
277  G4Exception("G4ParameterisationTrdX::CheckParametersValidity()",
278  "GeomDiv0001", FatalException, message);
279  }
280 */
281 }
282 
283 //--------------------------------------------------------------------------
286 {
287 }
288 
289 //--------------------------------------------------------------------------
292  G4double width, G4double offset,
293  G4VSolid* msolid, DivisionType divType)
294  : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
295 {
297  SetType( "DivisionTrdY" );
298 
299  G4Trd* msol = (G4Trd*)(fmotherSolid);
300  if( divType == DivWIDTH )
301  {
302  fnDiv = CalculateNDiv( 2*msol->GetYHalfLength1(),
303  width, offset );
304  }
305  else if( divType == DivNDIV )
306  {
308  nDiv, offset );
309  }
310 
311 #ifdef G4DIVDEBUG
312  if( verbose >= 1 )
313  {
314  G4cout << " G4ParameterisationTrdY no divisions " << fnDiv
315  << " = " << nDiv << G4endl
316  << " Offset " << foffset << " = " << offset << G4endl
317  << " width " << fwidth << " = " << width << G4endl;
318  }
319 #endif
320 }
321 
322 //------------------------------------------------------------------------
324 {
325 }
326 
327 //------------------------------------------------------------------------
329 {
330  G4Trd* msol = (G4Trd*)(fmotherSolid);
331  return 2*msol->GetYHalfLength1();
332 }
333 
334 //--------------------------------------------------------------------------
335 void
337 ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol ) const
338 {
339  G4Trd* msol = (G4Trd*)(fmotherSolid );
340  G4double mdy = msol->GetYHalfLength1();
341 
342  //----- translation
343  G4ThreeVector origin(0.,0.,0.);
344  G4double posi = -mdy + foffset + (copyNo+0.5)*fwidth;
345 
346  if( faxis == kYAxis )
347  {
348  origin.setY( posi );
349  }
350  else
351  {
352  std::ostringstream message;
353  message << "Only axes along Y are allowed ! Axis: " << faxis;
354  G4Exception("G4ParameterisationTrdY::ComputeTransformation()",
355  "GeomDiv0002", FatalException, message);
356  }
357 
358 #ifdef G4DIVDEBUG
359  if( verbose >= 2 )
360  {
361  G4cout << std::setprecision(8)
362  << " G4ParameterisationTrdY::ComputeTransformation " << copyNo
363  << " pos " << origin << " rot mat " << " axis " << faxis << G4endl;
364  }
365 #endif
366 
367  //----- set translation
368  physVol->SetTranslation( origin );
369 }
370 
371 //--------------------------------------------------------------------------
372 void
374 ComputeDimensions(G4Trd& trd, const G4int, const G4VPhysicalVolume*) const
375 {
376  //---- The division along Y of a Trd will result a Trd, only
377  //--- if Y at -Z and +Z are equal, else use the G4Trap version
378  G4Trd* msol = (G4Trd*)(fmotherSolid);
379 
380  G4double pDx1 = msol->GetXHalfLength1();
381  G4double pDx2 = msol->GetXHalfLength2();
382  G4double pDz = msol->GetZHalfLength();
383  G4double pDy = fwidth/2. - fhgap;
384 
385  trd.SetAllParameters ( pDx1, pDx2, pDy, pDy, pDz );
386 
387 #ifdef G4DIVDEBUG
388  if( verbose >= 2 )
389  {
390  G4cout << " G4ParameterisationTrdY::ComputeDimensions():" << G4endl;
391  trd.DumpInfo();
392  }
393 #endif
394 }
395 
396 //--------------------------------------------------------------------------
398 {
400 
401  G4Trd* msol = (G4Trd*)(fmotherSolid);
402 
403  G4double mpDy1 = msol->GetYHalfLength1();
404  G4double mpDy2 = msol->GetYHalfLength2();
405 
406  if( std::fabs(mpDy1 - mpDy2) > kCarTolerance )
407  {
408  std::ostringstream message;
409  message << "Invalid solid specification. NOT supported." << G4endl
410  << "Making a division of a TRD along axis Y while" << G4endl
411  << "the Y half lengths are not equal is not (yet)" << G4endl
412  << "supported. It will result in non-equal" << G4endl
413  << "division solids.";
414  G4Exception("G4ParameterisationTrdY::CheckParametersValidity()",
415  "GeomDiv0001", FatalException, message);
416  }
417 }
418 
419 //--------------------------------------------------------------------------
422  G4double width, G4double offset,
423  G4VSolid* msolid, DivisionType divType )
424  : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
425 {
427  SetType( "DivTrdZ" );
428 
429  G4Trd* msol = (G4Trd*)(fmotherSolid);
430  if( divType == DivWIDTH )
431  {
432  fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(),
433  width, offset );
434  }
435  else if( divType == DivNDIV )
436  {
437  fwidth = CalculateWidth( 2*msol->GetZHalfLength(),
438  nDiv, offset );
439  }
440 
441 #ifdef G4DIVDEBUG
442  if( verbose >= 1 )
443  {
444  G4cout << " G4ParameterisationTrdZ no divisions " << fnDiv
445  << " = " << nDiv << G4endl
446  << " Offset " << foffset << " = " << offset << G4endl
447  << " Width " << fwidth << " = " << width << G4endl;
448  }
449 #endif
450 }
451 
452 //------------------------------------------------------------------------
454 {
455 }
456 
457 //------------------------------------------------------------------------
459 {
460  G4Trd* msol = (G4Trd*)(fmotherSolid);
461  return 2*msol->GetZHalfLength();
462 }
463 
464 //--------------------------------------------------------------------------
465 void
467 ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
468 {
469  G4Trd* msol = (G4Trd*)(fmotherSolid );
470  G4double mdz = msol->GetZHalfLength();
471 
472  //----- translation
473  G4ThreeVector origin(0.,0.,0.);
474  G4double posi = -mdz + OffsetZ() + (copyNo+0.5)*fwidth;
475  if( faxis == kZAxis )
476  {
477  origin.setZ( posi );
478  }
479  else
480  {
481  std::ostringstream message;
482  message << "Only axes along Z are allowed ! Axis: " << faxis;
483  G4Exception("G4ParameterisationTrdZ::ComputeTransformation()",
484  "GeomDiv0002", FatalException, message);
485  }
486 
487 #ifdef G4DIVDEBUG
488  if( verbose >= 1 )
489  {
490  G4cout << std::setprecision(8) << " G4ParameterisationTrdZ: "
491  << copyNo << G4endl
492  << " Position: " << origin << " - Offset: " << foffset
493  << " - Width: " << fwidth << " Axis " << faxis << G4endl;
494  }
495 #endif
496 
497  //----- set translation
498  physVol->SetTranslation( origin );
499 }
500 
501 //--------------------------------------------------------------------------
502 void
504 ComputeDimensions(G4Trd& trd, const G4int copyNo,
505  const G4VPhysicalVolume*) const
506 {
507  //---- The division along Z of a Trd will result a Trd
508  G4Trd* msol = (G4Trd*)(fmotherSolid);
509 
510  G4double pDx1 = msol->GetXHalfLength1();
511  G4double DDx = (msol->GetXHalfLength2() - msol->GetXHalfLength1() );
512  G4double pDy1 = msol->GetYHalfLength1();
513  G4double DDy = (msol->GetYHalfLength2() - msol->GetYHalfLength1() );
514  G4double pDz = fwidth/2. - fhgap;
515  G4double zLength = 2*msol->GetZHalfLength();
516 
517  trd.SetAllParameters( pDx1+DDx*(OffsetZ()+copyNo*fwidth+fhgap)/zLength,
518  pDx1+DDx*(OffsetZ()+(copyNo+1)*fwidth-fhgap)/zLength,
519  pDy1+DDy*(OffsetZ()+copyNo*fwidth+fhgap)/zLength,
520  pDy1+DDy*(OffsetZ()+(copyNo+1)*fwidth-fhgap)/zLength,
521  pDz );
522 
523 #ifdef G4DIVDEBUG
524  if( verbose >= 1 )
525  {
526  G4cout << " G4ParameterisationTrdZ::ComputeDimensions()"
527  << " - Mother TRD " << G4endl;
528  msol->DumpInfo();
529  G4cout << " - Parameterised TRD: "
530  << copyNo << G4endl;
531  trd.DumpInfo();
532  }
533 #endif
534 }
void SetAllParameters(G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition: G4Trd.cc:168
G4String GetName() const
void ComputeDimensions(G4Trd &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const
G4double GetYHalfLength1() const
G4ParameterisationTrdX(EAxis axis, G4int nCopies, G4double width, G4double offset, G4VSolid *motherSolid, DivisionType divType)
CLHEP::Hep3Vector G4ThreeVector
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
G4ParameterisationTrdY(EAxis axis, G4int nCopies, G4double width, G4double offset, G4VSolid *motherSolid, DivisionType divType)
void SetType(const G4String &type)
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
#define width
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
Definition: G4Trd.hh:72
virtual G4GeometryType GetEntityType() const =0
G4double GetZHalfLength() const
int G4int
Definition: G4Types.hh:78
void DumpInfo() const
void ComputeDimensions(G4Trd &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const
G4double GetXHalfLength2() const
G4double GetMaxParameter() const
G4GLOB_DLL std::ostream G4cout
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
void SetAllParameters(G4double pDz, G4double pTheta, G4double pPhi, G4double pDy1, G4double pDx1, G4double pDx2, G4double pAlp1, G4double pDy2, G4double pDx3, G4double pDx4, G4double pAlp2)
Definition: G4Trap.cc:604
void SetTranslation(const G4ThreeVector &v)
G4double GetYHalfLength2() const
G4VParameterisationTrd(EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *msolid, DivisionType divType)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
EAxis
Definition: geomdefs.hh:54
G4double GetMaxParameter() const
G4double GetMaxParameter() const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
#define G4endl
Definition: G4ios.hh:61
G4double GetXHalfLength1() const
double G4double
Definition: G4Types.hh:76
G4ParameterisationTrdZ(EAxis axis, G4int nCopies, G4double width, G4double offset, G4VSolid *motherSolid, DivisionType divType)
void ComputeDimensions(G4Trd &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const