Geant4-11
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// G4ParameterisationTrd[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#include "G4ThreeVector.hh"
37#include "G4RotationMatrix.hh"
38#include "G4VPhysicalVolume.hh"
39#include "G4LogicalVolume.hh"
40#include "G4ReflectedSolid.hh"
41#include "G4Trd.hh"
42#include "G4Trap.hh"
43
44//--------------------------------------------------------------------------
47 G4double offset, G4VSolid* msolid,
48 DivisionType divType )
49 : G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid )
50{
51 G4Trd* msol = (G4Trd*)(msolid);
52 if (msolid->GetEntityType() == "G4ReflectedSolid")
53 {
54 // Get constituent solid
55 G4VSolid* mConstituentSolid
56 = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
57 msol = (G4Trd*)(mConstituentSolid);
58
59 // Create a new solid with inversed parameters
60 G4Trd* newSolid
61 = new G4Trd(msol->GetName(),
62 msol->GetXHalfLength2(), msol->GetXHalfLength1(),
63 msol->GetYHalfLength2(), msol->GetYHalfLength1(),
64 msol->GetZHalfLength());
65 msol = newSolid;
66 fmotherSolid = newSolid;
67 fReflectedSolid = true;
68 fDeleteSolid = true;
69 }
70}
71
72//------------------------------------------------------------------------
74{
75}
76
77//------------------------------------------------------------------------
80 G4double width, G4double offset,
81 G4VSolid* msolid, DivisionType divType )
82 : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
83{
85 SetType( "DivisionTrdX" );
86
87 G4Trd* msol = (G4Trd*)(fmotherSolid);
88 if( divType == DivWIDTH )
89 {
91 width, offset );
92 }
93 else if( divType == DivNDIV )
94 {
96 nDiv, offset );
97 }
98
99#ifdef G4DIVDEBUG
100 if( verbose >= 1 )
101 {
102 G4cout << " G4ParameterisationTrdX - ## divisions " << fnDiv << " = "
103 << nDiv << G4endl
104 << " Offset " << foffset << " = " << offset << G4endl
105 << " Width " << fwidth << " = " << width << G4endl;
106 }
107#endif
108
109 G4double mpDx1 = msol->GetXHalfLength1();
110 G4double mpDx2 = msol->GetXHalfLength2();
111 if( std::fabs(mpDx1 - mpDx2) > kCarTolerance )
112 {
113 bDivInTrap = true;
114 }
115}
116
117//------------------------------------------------------------------------
119{
120}
121
122//------------------------------------------------------------------------
124{
125 G4Trd* msol = (G4Trd*)(fmotherSolid);
126 return (msol->GetXHalfLength1()+msol->GetXHalfLength2());
127}
128
129//------------------------------------------------------------------------
130void
132ComputeTransformation( const G4int copyNo,
133 G4VPhysicalVolume *physVol ) const
134{
135 G4Trd* msol = (G4Trd*)(fmotherSolid );
136 G4double mdx = ( msol->GetXHalfLength1() + msol->GetXHalfLength2() ) / 2.;
137 //----- translation
138 G4ThreeVector origin(0.,0.,0.);
139 G4double posi;
140 posi = -mdx + foffset + (copyNo+0.5)*fwidth;
141 if( faxis == kXAxis )
142 {
143 origin.setX( posi );
144 }
145 else
146 {
147 std::ostringstream message;
148 message << "Only axes along X are allowed ! Axis: " << faxis;
149 G4Exception("G4ParameterisationTrdX::ComputeTransformation()",
150 "GeomDiv0002", FatalException, message);
151 }
152
153#ifdef G4DIVDEBUG
154 if( verbose >= 2 )
155 {
156 G4cout << std::setprecision(8)
157 << " G4ParameterisationTrdX::ComputeTransformation() "
158 << copyNo << G4endl
159 << " Position: " << origin << " - Axis: " << faxis << G4endl;
160 }
161#endif
162
163 //----- set translation
164 physVol->SetTranslation( origin );
165}
166
167//--------------------------------------------------------------------------
168void
170ComputeDimensions( G4Trd& trd, [[maybe_unused]] const G4int copyNo, const G4VPhysicalVolume* ) const
171{
172 G4Trd* msol = (G4Trd*)(fmotherSolid);
173 G4double pDy1 = msol->GetYHalfLength1();
174 G4double pDy2 = msol->GetYHalfLength2();
175 G4double pDz = msol->GetZHalfLength();
176 G4double pDx = fwidth/2. - fhgap;
177
178 trd.SetAllParameters ( pDx, pDx, pDy1, pDy2, pDz );
179
180#ifdef G4DIVDEBUG
181 if( verbose >= 2 )
182 {
183 G4cout << " G4ParameterisationTrdX::ComputeDimensions():"
184 << copyNo << G4endl;
185 trd.DumpInfo();
186 }
187#endif
188}
189
190//--------------------------------------------------------------------------
191void
193 const G4VPhysicalVolume* ) const
194{
195 G4Trd* msol = (G4Trd*)(fmotherSolid);
196 G4double pDy1 = msol->GetYHalfLength1();
197 G4double pDy2 = msol->GetYHalfLength2();
198 G4double pDz = msol->GetZHalfLength();
199 //fwidth is at Z=0
200 G4double pDx1 = msol->GetXHalfLength1();
201 G4double pDx2 = msol->GetXHalfLength2();
202 // G4double pDxAVE = (pDx1+pDx2)/2.;
203 G4double xChangeRatio = (pDx2-pDx1) / (pDx2+pDx1);
204 G4double fWidChange = xChangeRatio*fwidth;
205 G4double fWid1 = fwidth - fWidChange;
206 G4double fWid2 = fwidth + fWidChange;
207 G4double fOffsetChange = xChangeRatio*foffset/2.;
208 G4double fOffset1 = foffset - fOffsetChange;
209 G4double fOffset2 = foffset + fOffsetChange;
210 G4double cxy1 = -pDx1+fOffset1 + (copyNo+0.5)*fWid1;// centre of the side at y=-pDy1
211 G4double cxy2 = -pDx2+fOffset2 + (copyNo+0.5)*fWid2;// centre of the side at y=+pDy1
212 G4double alp = std::atan( (cxy2-cxy1)/(pDz*2.) );
213
214 pDx1 = fwidth/2. - fWidChange/2.;
215 pDx2 = fwidth/2. + fWidChange/2.;
216
217
218 trap.SetAllParameters ( pDz,
219 alp,
220 0.,
221 pDy1,
222 pDx1,
223 pDx1,
224 0.,
225 pDy2,
226 pDx2 - fhgap,
227 pDx2 - fhgap * pDx2/pDx1,
228 0.);
229
230#ifdef G4DIVDEBUG
231 if( verbose >= 2 )
232 {
233 G4cout << " G4ParameterisationTrdX::ComputeDimensions():"
234 << copyNo << G4endl;
235 trap.DumpInfo();
236 }
237#endif
238}
239
240
241//--------------------------------------------------------------------------
244 G4double width, G4double offset,
245 G4VSolid* msolid, DivisionType divType)
246 : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
247{
249 SetType( "DivisionTrdY" );
250
251 G4Trd* msol = (G4Trd*)(fmotherSolid);
252 if( divType == DivWIDTH )
253 {
255 width, offset );
256 }
257 else if( divType == DivNDIV )
258 {
260 nDiv, offset );
261 }
262
263#ifdef G4DIVDEBUG
264 if( verbose >= 1 )
265 {
266 G4cout << " G4ParameterisationTrdY no divisions " << fnDiv
267 << " = " << nDiv << G4endl
268 << " Offset " << foffset << " = " << offset << G4endl
269 << " width " << fwidth << " = " << width << G4endl;
270 }
271#endif
272}
273
274//------------------------------------------------------------------------
276{
277}
278
279//------------------------------------------------------------------------
281{
282 G4Trd* msol = (G4Trd*)(fmotherSolid);
283 return (msol->GetYHalfLength1()+msol->GetYHalfLength2());
284}
285
286//--------------------------------------------------------------------------
287void
289ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol ) const
290{
291 G4Trd* msol = (G4Trd*)(fmotherSolid );
292 G4double mdy = ( msol->GetYHalfLength1() + msol->GetYHalfLength2() ) / 2.;
293
294 //----- translation
295 G4ThreeVector origin(0.,0.,0.);
296 G4double posi = -mdy + foffset + (copyNo+0.5)*fwidth;
297
298 if( faxis == kYAxis )
299 {
300 origin.setY( posi );
301 }
302 else
303 {
304 std::ostringstream message;
305 message << "Only axes along Y are allowed ! Axis: " << faxis;
306 G4Exception("G4ParameterisationTrdY::ComputeTransformation()",
307 "GeomDiv0002", FatalException, message);
308 }
309
310#ifdef G4DIVDEBUG
311 if( verbose >= 2 )
312 {
313 G4cout << std::setprecision(8)
314 << " G4ParameterisationTrdY::ComputeTransformation " << copyNo
315 << " pos " << origin << " rot mat " << " axis " << faxis << G4endl;
316 }
317#endif
318
319 //----- set translation
320 physVol->SetTranslation( origin );
321}
322
323//--------------------------------------------------------------------------
324void
326ComputeDimensions(G4Trd& trd, const G4int, const G4VPhysicalVolume*) const
327{
328 //---- The division along Y of a Trd will result a Trd, only
329 //--- if Y at -Z and +Z are equal, else use the G4Trap version
330 G4Trd* msol = (G4Trd*)(fmotherSolid);
331
332 G4double pDx1 = msol->GetXHalfLength1();
333 G4double pDx2 = msol->GetXHalfLength2();
334 G4double pDz = msol->GetZHalfLength();
335 G4double pDy = fwidth/2. - fhgap;
336
337 trd.SetAllParameters ( pDx1, pDx2, pDy, pDy, pDz );
338
339#ifdef G4DIVDEBUG
340 if( verbose >= 2 )
341 {
342 G4cout << " G4ParameterisationTrdY::ComputeDimensions():" << G4endl;
343 trd.DumpInfo();
344 }
345#endif
346}
347
348//--------------------------------------------------------------------------
349void
351 const G4VPhysicalVolume* ) const
352{
353 G4Trd* msol = (G4Trd*)(fmotherSolid);
354 G4double pDx1 = msol->GetXHalfLength1();
355 G4double pDx2 = msol->GetXHalfLength2();
356 G4double pDz = msol->GetZHalfLength();
357 //fwidth is at Z=0
358 G4double pDy1 = msol->GetYHalfLength1();
359 G4double pDy2 = msol->GetYHalfLength2();
360 // G4double pDxAVE = (pDy1+pDy2)/2.;
361 G4double yChangeRatio = (pDy2-pDy1) / (pDy2+pDy1);
362 G4double fWidChange = yChangeRatio*fwidth;
363 G4double fWid1 = fwidth - fWidChange;
364 G4double fWid2 = fwidth + fWidChange;
365 G4double fOffsetChange = yChangeRatio*foffset/2.;
366 G4double fOffset1 = foffset - fOffsetChange;
367 G4double fOffset2 = foffset + fOffsetChange;
368 G4double cyx1 = -pDy1+fOffset1 + (copyNo+0.5)*fWid1;// centre of the side at y=-pDy1
369 G4double cyx2 = -pDy2+fOffset2 + (copyNo+0.5)*fWid2;// centre of the side at y=+pDy1
370 G4double alp = std::atan( (cyx2-cyx1)/(pDz*2.) );
371
372 pDy1 = fwidth/2. - fWidChange/2.;
373 pDy2 = fwidth/2. + fWidChange/2.;
374
375
376 trap.SetAllParameters ( pDz,
377 alp,
378 90*CLHEP::deg,
379 pDy1,
380 pDx1,
381 pDx1,
382 0.,
383 pDy2,
384 pDx2 - fhgap,
385 pDx2 - fhgap * pDx2/pDx1,
386 0.);
387
388#ifdef G4DIVDEBUG
389 if( verbose >= 2 )
390 {
391 G4cout << " G4ParameterisationTrdY::ComputeDimensions():"
392 << copyNo << G4endl;
393 trap.DumpInfo();
394 }
395#endif
396}
397
398//--------------------------------------------------------------------------
401 G4double width, G4double offset,
402 G4VSolid* msolid, DivisionType divType )
403 : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType )
404{
406 SetType( "DivTrdZ" );
407
408 G4Trd* msol = (G4Trd*)(fmotherSolid);
409 if( divType == DivWIDTH )
410 {
412 width, offset );
413 }
414 else if( divType == DivNDIV )
415 {
417 nDiv, offset );
418 }
419
420#ifdef G4DIVDEBUG
421 if( verbose >= 1 )
422 {
423 G4cout << " G4ParameterisationTrdZ no divisions " << fnDiv
424 << " = " << nDiv << G4endl
425 << " Offset " << foffset << " = " << offset << G4endl
426 << " Width " << fwidth << " = " << width << G4endl;
427 }
428#endif
429}
430
431//------------------------------------------------------------------------
433{
434}
435
436//------------------------------------------------------------------------
438{
439 G4Trd* msol = (G4Trd*)(fmotherSolid);
440 return 2*msol->GetZHalfLength();
441}
442
443//--------------------------------------------------------------------------
444void
446ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
447{
448 G4Trd* msol = (G4Trd*)(fmotherSolid );
449 G4double mdz = msol->GetZHalfLength();
450
451 //----- translation
452 G4ThreeVector origin(0.,0.,0.);
453 G4double posi = -mdz + OffsetZ() + (copyNo+0.5)*fwidth;
454 if( faxis == kZAxis )
455 {
456 origin.setZ( posi );
457 }
458 else
459 {
460 std::ostringstream message;
461 message << "Only axes along Z are allowed ! Axis: " << faxis;
462 G4Exception("G4ParameterisationTrdZ::ComputeTransformation()",
463 "GeomDiv0002", FatalException, message);
464 }
465
466#ifdef G4DIVDEBUG
467 if( verbose >= 1 )
468 {
469 G4cout << std::setprecision(8) << " G4ParameterisationTrdZ: "
470 << copyNo << G4endl
471 << " Position: " << origin << " - Offset: " << foffset
472 << " - Width: " << fwidth << " Axis " << faxis << G4endl;
473 }
474#endif
475
476 //----- set translation
477 physVol->SetTranslation( origin );
478}
479
480//--------------------------------------------------------------------------
481void
483ComputeDimensions(G4Trd& trd, const G4int copyNo,
484 const G4VPhysicalVolume*) const
485{
486 //---- The division along Z of a Trd will result a Trd
487 G4Trd* msol = (G4Trd*)(fmotherSolid);
488
489 G4double pDx1 = msol->GetXHalfLength1();
490 G4double DDx = (msol->GetXHalfLength2() - msol->GetXHalfLength1() );
491 G4double pDy1 = msol->GetYHalfLength1();
492 G4double DDy = (msol->GetYHalfLength2() - msol->GetYHalfLength1() );
493 G4double pDz = fwidth/2. - fhgap;
494 G4double zLength = 2*msol->GetZHalfLength();
495
496 trd.SetAllParameters( pDx1+DDx*(OffsetZ()+copyNo*fwidth+fhgap)/zLength,
497 pDx1+DDx*(OffsetZ()+(copyNo+1)*fwidth-fhgap)/zLength,
498 pDy1+DDy*(OffsetZ()+copyNo*fwidth+fhgap)/zLength,
499 pDy1+DDy*(OffsetZ()+(copyNo+1)*fwidth-fhgap)/zLength,
500 pDz );
501
502#ifdef G4DIVDEBUG
503 if( verbose >= 1 )
504 {
505 G4cout << " G4ParameterisationTrdZ::ComputeDimensions()"
506 << " - Mother TRD " << G4endl;
507 msol->DumpInfo();
508 G4cout << " - Parameterised TRD: "
509 << copyNo << G4endl;
510 trd.DumpInfo();
511 }
512#endif
513}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void setY(double)
void setZ(double)
void setX(double)
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
void ComputeDimensions(G4Trd &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const
G4ParameterisationTrdX(EAxis axis, G4int nCopies, G4double width, G4double offset, G4VSolid *motherSolid, DivisionType divType)
G4ParameterisationTrdY(EAxis axis, G4int nCopies, G4double width, G4double offset, G4VSolid *motherSolid, DivisionType divType)
void ComputeDimensions(G4Trd &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
void ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const
void ComputeDimensions(G4Trd &trd, const G4int copyNo, const G4VPhysicalVolume *pv) const
G4ParameterisationTrdZ(EAxis axis, G4int nCopies, G4double width, G4double offset, G4VSolid *motherSolid, DivisionType divType)
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:282
Definition: G4Trd.hh:63
void SetAllParameters(G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition: G4Trd.cc:128
G4double GetXHalfLength2() const
G4double GetYHalfLength2() const
G4double GetXHalfLength1() const
G4double GetYHalfLength1() const
G4double GetZHalfLength() const
void SetType(const G4String &type)
G4double CalculateWidth(G4double motherDim, G4int nDiv, G4double offset) const
G4int CalculateNDiv(G4double motherDim, G4double width, G4double offset) const
G4VParameterisationTrd(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
@ kYAxis
Definition: geomdefs.hh:56
@ kXAxis
Definition: geomdefs.hh:55
@ kZAxis
Definition: geomdefs.hh:57
static constexpr double deg