source: trunk/source/geometry/divisions/src/G4ParameterisationPolyhedra.cc@ 1296

Last change on this file since 1296 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 25.8 KB
Line 
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: G4ParameterisationPolyhedra.cc,v 1.19 2009/05/20 08:35:52 ivana Exp $
28// GEANT4 tag $Name: geant4-09-03 $
29//
30// class G4ParameterisationPolyhedra Implementation file
31//
32// 14.10.03 - P.Arce, Initial version
33// 08.04.04 - I.Hrivnacova, Implemented reflection
34// --------------------------------------------------------------------
35
36#include "G4ParameterisationPolyhedra.hh"
37
38#include <iomanip>
39#include "G4ThreeVector.hh"
40#include "G4GeometryTolerance.hh"
41#include "G4RotationMatrix.hh"
42#include "G4VPhysicalVolume.hh"
43#include "G4LogicalVolume.hh"
44#include "G4ReflectedSolid.hh"
45#include "G4Polyhedra.hh"
46
47//--------------------------------------------------------------------------
48G4VParameterisationPolyhedra::
49G4VParameterisationPolyhedra( EAxis axis, G4int nDiv, G4double width,
50 G4double offset, G4VSolid* msolid,
51 DivisionType divType )
52 : G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid )
53{
54 G4Polyhedra* msol = (G4Polyhedra*)(msolid);
55 if ((msolid->GetEntityType() != "G4ReflectedSolid") && (msol->IsGeneric()))
56 {
57 G4String message =
58 "Sorry, generic construct for G4Polyhedra NOT supported.\n Solid: "
59 + msol->GetName();
60 G4Exception("G4VParameterisationPolyhedra::G4VParameterisationPolyhedra()",
61 "NotSupported", FatalException, message);
62 }
63 if (msolid->GetEntityType() == "G4ReflectedSolid")
64 {
65 // Get constituent solid
66 G4VSolid* mConstituentSolid
67 = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
68 msol = (G4Polyhedra*)(mConstituentSolid);
69
70 // Get parameters
71 G4int nofSides = msol->GetOriginalParameters()->numSide;
72 G4int nofZplanes = msol->GetOriginalParameters()->Num_z_planes;
73 G4double* zValues = msol->GetOriginalParameters()->Z_values;
74 G4double* rminValues = msol->GetOriginalParameters()->Rmin;
75 G4double* rmaxValues = msol->GetOriginalParameters()->Rmax;
76
77 // Invert z values,
78 // convert radius parameters
79 G4double* rminValues2 = new G4double[nofZplanes];
80 G4double* rmaxValues2 = new G4double[nofZplanes];
81 G4double* zValuesRefl = new G4double[nofZplanes];
82 for (G4int i=0; i<nofZplanes; i++)
83 {
84 rminValues2[i] = rminValues[i] * ConvertRadiusFactor(*msol);
85 rmaxValues2[i] = rmaxValues[i] * ConvertRadiusFactor(*msol);
86 zValuesRefl[i] = - zValues[i];
87 }
88
89 G4Polyhedra* newSolid
90 = new G4Polyhedra(msol->GetName(),
91 msol->GetStartPhi(),
92 msol->GetEndPhi() - msol->GetStartPhi(),
93 nofSides,
94 nofZplanes, zValuesRefl, rminValues2, rmaxValues2);
95
96 delete [] rminValues2;
97 delete [] rmaxValues2;
98 delete [] zValuesRefl;
99
100 msol = newSolid;
101 fmotherSolid = newSolid;
102 fReflectedSolid = true;
103 fDeleteSolid = true;
104 }
105}
106
107//------------------------------------------------------------------------
108G4VParameterisationPolyhedra::~G4VParameterisationPolyhedra()
109{
110}
111
112//--------------------------------------------------------------------------
113G4double
114G4VParameterisationPolyhedra::
115ConvertRadiusFactor(const G4Polyhedra& phedra) const
116{
117 G4double phiTotal = phedra.GetEndPhi() - phedra.GetStartPhi();
118 G4int nofSides = phedra.GetOriginalParameters()->numSide;
119
120 if ( (phiTotal <=0) || (phiTotal >
121 2*pi+G4GeometryTolerance::GetInstance()->GetAngularTolerance()) )
122 { phiTotal = 2*pi; }
123
124 return std::cos(0.5*phiTotal/nofSides);
125}
126
127//--------------------------------------------------------------------------
128G4ParameterisationPolyhedraRho::
129G4ParameterisationPolyhedraRho( EAxis axis, G4int nDiv,
130 G4double width, G4double offset,
131 G4VSolid* msolid, DivisionType divType )
132 : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType )
133{
134
135 CheckParametersValidity();
136 SetType( "DivisionPolyhedraRho" );
137
138 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
139 G4PolyhedraHistorical* original_pars = msol->GetOriginalParameters();
140
141 if( divType == DivWIDTH )
142 {
143 fnDiv = CalculateNDiv( original_pars->Rmax[0]
144 - original_pars->Rmin[0], width, offset );
145 }
146 else if( divType == DivNDIV )
147 {
148 fwidth = CalculateWidth( original_pars->Rmax[0]
149 - original_pars->Rmin[0], nDiv, offset );
150 }
151
152#ifdef G4DIVDEBUG
153 if( verbose >= 1 )
154 {
155 G4cout << " G4ParameterisationPolyhedraRho - # divisions " << fnDiv
156 << " = " << nDiv << G4endl
157 << " Offset " << foffset << " = " << offset << G4endl
158 << " Width " << fwidth << " = " << width << G4endl;
159 }
160#endif
161}
162
163//------------------------------------------------------------------------
164G4ParameterisationPolyhedraRho::~G4ParameterisationPolyhedraRho()
165{
166}
167
168//---------------------------------------------------------------------
169void G4ParameterisationPolyhedraRho::CheckParametersValidity()
170{
171 G4VDivisionParameterisation::CheckParametersValidity();
172
173 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
174
175 if( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH )
176 {
177 G4cerr << "WARNING - "
178 << "G4ParameterisationPolyhedraRho::CheckParametersValidity()"
179 << G4endl
180 << " Solid " << msol->GetName() << G4endl
181 << " Division along R will be done with a width "
182 << "different for each solid section." << G4endl
183 << " WIDTH will not be used !" << G4endl;
184 }
185 if( foffset != 0. )
186 {
187 G4cerr << "WARNING - "
188 << "G4ParameterisationPolyhedraRho::CheckParametersValidity()"
189 << G4endl
190 << " Solid " << msol->GetName() << G4endl
191 << " Division along R will be done with a width "
192 << "different for each solid section." << G4endl
193 << " OFFSET will not be used !" << G4endl;
194 }
195}
196
197//------------------------------------------------------------------------
198G4double G4ParameterisationPolyhedraRho::GetMaxParameter() const
199{
200 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
201 G4PolyhedraHistorical* original_pars = msol->GetOriginalParameters();
202 return original_pars->Rmax[0] - original_pars->Rmin[0];
203}
204
205//--------------------------------------------------------------------------
206void
207G4ParameterisationPolyhedraRho::
208ComputeTransformation( const G4int, G4VPhysicalVolume* physVol ) const
209{
210 //----- translation
211 G4ThreeVector origin(0.,0.,0.);
212
213 //----- set translation
214 physVol->SetTranslation( origin );
215
216 //----- calculate rotation matrix: unit
217
218#ifdef G4DIVDEBUG
219 if( verbose >= 2 )
220 {
221 G4cout << " G4ParameterisationPolyhedraRho " << G4endl
222 << " foffset: " << foffset/deg
223 << " - fwidth: " << fwidth/deg << G4endl;
224 }
225#endif
226
227 ChangeRotMatrix( physVol );
228
229#ifdef G4DIVDEBUG
230 if( verbose >= 2 )
231 {
232 G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraRho "
233 << G4endl
234 << " Position: " << origin
235 << " - Width: " << fwidth
236 << " - Axis: " << faxis << G4endl;
237 }
238#endif
239}
240
241//--------------------------------------------------------------------------
242void
243G4ParameterisationPolyhedraRho::
244ComputeDimensions( G4Polyhedra& phedra, const G4int copyNo,
245 const G4VPhysicalVolume* ) const
246{
247 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
248
249 G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters();
250 G4PolyhedraHistorical origparam( *origparamMother );
251 G4int nZplanes = origparamMother->Num_z_planes;
252
253 G4double width = 0.;
254 for( G4int ii = 0; ii < nZplanes; ii++ )
255 {
256 width = CalculateWidth( origparamMother->Rmax[ii]
257 - origparamMother->Rmin[ii], fnDiv, foffset );
258 origparam.Rmin[ii] = origparamMother->Rmin[ii]+foffset+width*copyNo;
259 origparam.Rmax[ii] = origparamMother->Rmin[ii]+foffset+width*(copyNo+1);
260 }
261
262 phedra.SetOriginalParameters(&origparam); // copy values & transfer pointers
263 phedra.Reset(); // reset to new solid parameters
264
265#ifdef G4DIVDEBUG
266 if( verbose >= -2 )
267 {
268 G4cout << "G4ParameterisationPolyhedraRho::ComputeDimensions()" << G4endl
269 << "-- Parametrised phedra copy-number: " << copyNo << G4endl;
270 phedra.DumpInfo();
271 }
272#endif
273}
274
275//--------------------------------------------------------------------------
276G4ParameterisationPolyhedraPhi::
277G4ParameterisationPolyhedraPhi( EAxis axis, G4int nDiv,
278 G4double width, G4double offset,
279 G4VSolid* msolid, DivisionType divType )
280 : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType )
281{
282 CheckParametersValidity();
283 SetType( "DivisionPolyhedraPhi" );
284
285 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
286 G4double deltaPhi = msol->GetEndPhi() - msol->GetStartPhi();
287
288 if( divType == DivWIDTH )
289 {
290 fnDiv = msol->GetNumSide();
291 }
292
293 fwidth = CalculateWidth( deltaPhi, fnDiv, 0.0 );
294
295#ifdef G4DIVDEBUG
296 if( verbose >= 1 )
297 {
298 G4cout << " G4ParameterisationPolyhedraPhi - # divisions " << fnDiv
299 << " = " << nDiv << G4endl
300 << " Offset " << foffset << " = " << offset << G4endl
301 << " Width " << fwidth << " = " << width << G4endl;
302 }
303#endif
304}
305
306//------------------------------------------------------------------------
307G4ParameterisationPolyhedraPhi::~G4ParameterisationPolyhedraPhi()
308{
309}
310
311//------------------------------------------------------------------------
312G4double G4ParameterisationPolyhedraPhi::GetMaxParameter() const
313{
314 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
315 return msol->GetEndPhi() - msol->GetStartPhi();
316}
317
318//---------------------------------------------------------------------
319void G4ParameterisationPolyhedraPhi::CheckParametersValidity()
320{
321 G4VDivisionParameterisation::CheckParametersValidity();
322
323 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
324
325 if( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH )
326 {
327 G4cerr << "WARNING - "
328 << "G4ParameterisationPolyhedraPhi::CheckParametersValidity()"
329 << G4endl
330 << " Solid " << msol->GetName() << G4endl
331 << " Division along PHI will be done splitting "
332 << "in the defined numSide." << G4endl
333 << " WIDTH will not be used !" << G4endl;
334 }
335 if( foffset != 0. )
336 {
337 G4cerr << "WARNING - "
338 << "G4ParameterisationPolyhedraPhi::CheckParametersValidity()"
339 << G4endl
340 << " Solid " << msol->GetName() << G4endl
341 << " Division along PHI will be done splitting "
342 << "in the defined numSide." << G4endl
343 << " OFFSET will not be used !" << G4endl;
344 }
345
346 G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters();
347
348 if( origparamMother->numSide != fnDiv && fDivisionType != DivWIDTH)
349 {
350 G4cerr << "ERROR - "
351 << "G4ParameterisationPolyhedraPhi::CheckParametersValidity()"
352 << G4endl
353 << " Division along PHI will be done splitting in the defined"
354 << G4endl
355 << " numSide, i.e, the number of division would be :"
356 << " " << origparamMother->numSide
357 << " instead of " << fnDiv << " !"
358 << G4endl;
359 G4Exception("G4ParameterisationPolyhedraPhi::CheckParametersValidity()",
360 "IllegalConstruct", FatalException,
361 "Not supported configuration.");
362 }
363}
364
365//--------------------------------------------------------------------------
366void
367G4ParameterisationPolyhedraPhi::
368ComputeTransformation( const G4int copyNo, G4VPhysicalVolume *physVol ) const
369{
370 //----- translation
371 G4ThreeVector origin(0.,0.,0.);
372 //----- set translation
373 physVol->SetTranslation( origin );
374
375 //----- calculate rotation matrix (so that all volumes point to the centre)
376 G4double posi = copyNo*fwidth;
377
378#ifdef G4DIVDEBUG
379 if( verbose >= 2 )
380 {
381 G4cout << " G4ParameterisationPolyhedraPhi - position: " << posi/deg
382 << G4endl
383 << " copyNo: " << copyNo
384 << " - fwidth: " << fwidth/deg << G4endl;
385 }
386#endif
387
388 ChangeRotMatrix( physVol, -posi );
389
390#ifdef G4DIVDEBUG
391 if( verbose >= 2 )
392 {
393 G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraPhi " << copyNo
394 << G4endl
395 << " Position: " << origin << " - Width: " << fwidth
396 << " - Axis: " << faxis << G4endl;
397 }
398#endif
399}
400
401//--------------------------------------------------------------------------
402void
403G4ParameterisationPolyhedraPhi::
404ComputeDimensions( G4Polyhedra& phedra, const G4int,
405 const G4VPhysicalVolume* ) const
406{
407 G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
408
409 G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters();
410 G4PolyhedraHistorical origparam( *origparamMother );
411
412 origparam.numSide = 1;
413 origparam.Start_angle = origparamMother->Start_angle;
414 origparam.Opening_angle = fwidth;
415
416 phedra.SetOriginalParameters(&origparam); // copy values & transfer pointers
417 phedra.Reset(); // reset to new solid parameters
418
419#ifdef G4DIVDEBUG
420 if( verbose >= 2 )
421 {
422 G4cout << "G4ParameterisationPolyhedraPhi::ComputeDimensions():" << G4endl;
423 phedra.DumpInfo();
424 }
425#endif
426}
427
428//--------------------------------------------------------------------------
429G4ParameterisationPolyhedraZ::
430G4ParameterisationPolyhedraZ( EAxis axis, G4int nDiv,
431 G4double width, G4double offset,
432 G4VSolid* msolid, DivisionType divType )
433 : G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType ),
434 fNSegment(0),
435 fOrigParamMother(((G4Polyhedra*)fmotherSolid)->GetOriginalParameters())
436{
437 CheckParametersValidity();
438 SetType( "DivisionPolyhedraZ" );
439
440 if( divType == DivWIDTH )
441 {
442 fnDiv =
443 CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
444 - fOrigParamMother->Z_values[0] , width, offset );
445 }
446 else if( divType == DivNDIV )
447 {
448 fwidth =
449 CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
450 - fOrigParamMother->Z_values[0] , nDiv, offset );
451 }
452
453#ifdef G4DIVDEBUG
454 if( verbose >= 1 )
455 {
456 G4cout << " G4ParameterisationPolyhedraZ - # divisions " << fnDiv << " = "
457 << nDiv << G4endl
458 << " Offset " << foffset << " = " << offset << G4endl
459 << " Width " << fwidth << " = " << width << G4endl;
460 }
461#endif
462}
463
464//---------------------------------------------------------------------
465G4ParameterisationPolyhedraZ::~G4ParameterisationPolyhedraZ()
466{
467}
468
469//------------------------------------------------------------------------
470G4double G4ParameterisationPolyhedraZ::GetR(G4double z,
471 G4double z1, G4double r1,
472 G4double z2, G4double r2) const
473{
474 // Linear parameterisation:
475 // r = az + b
476 // a = (r1 - r2)/(z1-z2)
477 // b = r1 - a*z1
478
479 return (r1-r2)/(z1-z2)*z + ( r1 - (r1-r2)/(z1-z2)*z1 ) ;
480}
481
482//------------------------------------------------------------------------
483G4double G4ParameterisationPolyhedraZ::GetRmin(G4double z, G4int nseg) const
484{
485// Get Rmin in the given z position for the given polyhedra segment
486
487 return GetR(z,
488 fOrigParamMother->Z_values[nseg],
489 fOrigParamMother->Rmin[nseg],
490 fOrigParamMother->Z_values[nseg+1],
491 fOrigParamMother->Rmin[nseg+1]);
492}
493
494//------------------------------------------------------------------------
495G4double G4ParameterisationPolyhedraZ::GetRmax(G4double z, G4int nseg) const
496{
497// Get Rmax in the given z position for the given polyhedra segment
498
499 return GetR(z,
500 fOrigParamMother->Z_values[nseg],
501 fOrigParamMother->Rmax[nseg],
502 fOrigParamMother->Z_values[nseg+1],
503 fOrigParamMother->Rmax[nseg+1]);
504}
505
506//------------------------------------------------------------------------
507G4double G4ParameterisationPolyhedraZ::GetMaxParameter() const
508{
509 return std::abs (fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
510 -fOrigParamMother->Z_values[0]);
511}
512
513//---------------------------------------------------------------------
514void G4ParameterisationPolyhedraZ::CheckParametersValidity()
515{
516 G4VDivisionParameterisation::CheckParametersValidity();
517
518 // Division will be following the mother polyhedra segments
519 if( fDivisionType == DivNDIV ) {
520 if( fOrigParamMother->Num_z_planes-1 != fnDiv ) {
521 G4cerr << "ERROR - "
522 << "G4ParameterisationPolyhedraZ::CheckParametersValidity()"
523 << G4endl
524 << " Division along Z will be done splitting in the defined"
525 << G4endl
526 << " z_planes, i.e, the number of division would be :"
527 << " " << fOrigParamMother->Num_z_planes-1
528 << " instead of " << fnDiv << " !"
529 << G4endl;
530 G4Exception("G4ParameterisationPolyhedraZ::CheckParametersValidity()",
531 "IllegalConstruct", FatalException,
532 "Not supported configuration.");
533 }
534 }
535
536 // Division will be done within one polyhedra segment
537 // with applying given width and offset
538 if( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH ) {
539 // Check if divided region does not span over more
540 // than one z segment
541
542 G4int isegstart = -1; // number of the segment containing start position
543 G4int isegend = -1; // number of the segment containing end position
544
545 if ( ! fReflectedSolid ) {
546 // The start/end position of the divided region
547 G4double zstart
548 = fOrigParamMother->Z_values[0] + foffset;
549 G4double zend
550 = fOrigParamMother->Z_values[0] + foffset + fnDiv* fwidth;
551
552 G4int counter = 0;
553 while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 ) {
554 // first segment
555 if ( zstart >= fOrigParamMother->Z_values[counter] &&
556 zstart < fOrigParamMother->Z_values[counter+1] ) {
557 isegstart = counter;
558 }
559 // last segment
560 if ( zend > fOrigParamMother->Z_values[counter] &&
561 zend <= fOrigParamMother->Z_values[counter+1] ) {
562 isegend = counter;
563 }
564 ++counter;
565 }
566 }
567 else {
568 // The start/end position of the divided region
569 G4double zstart
570 = fOrigParamMother->Z_values[0] - foffset;
571 G4double zend
572 = fOrigParamMother->Z_values[0] - ( foffset + fnDiv* fwidth);
573
574 G4int counter = 0;
575 while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 ) {
576 // first segment
577 if ( zstart <= fOrigParamMother->Z_values[counter] &&
578 zstart > fOrigParamMother->Z_values[counter+1] ) {
579 isegstart = counter;
580 }
581 // last segment
582 if ( zend < fOrigParamMother->Z_values[counter] &&
583 zend >= fOrigParamMother->Z_values[counter+1] ) {
584 isegend = counter;
585 }
586 ++counter;
587 }
588 }
589
590 if ( isegstart != isegend ) {
591 G4cerr << "WARNING - "
592 << "G4ParameterisationPolyhedraZ::CheckParametersValidity()"
593 << G4endl
594 << " Division with user defined width." << G4endl
595 << " Solid " << fmotherSolid->GetName() << G4endl
596 << " Divided region is not between two z planes."
597 << G4endl;
598
599 G4Exception("G4ParameterisationPolyhedraZ::CheckParametersValidity()",
600 "IllegalConstruct", FatalException,
601 "Not supported configuration.");
602 }
603
604 fNSegment = isegstart;
605 }
606}
607
608//---------------------------------------------------------------------
609void
610G4ParameterisationPolyhedraZ::
611ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol) const
612{
613 if ( fDivisionType == DivNDIV ) {
614 // The position of the centre of copyNo-th mother polycone segment
615 G4double posi = ( fOrigParamMother->Z_values[copyNo]
616 + fOrigParamMother->Z_values[copyNo+1])/2;
617 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
618 }
619
620 if ( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH ) {
621 // The position of the centre of copyNo-th division
622
623 G4double posi = fOrigParamMother->Z_values[0];
624
625 if ( ! fReflectedSolid )
626 posi += foffset + (2*copyNo + 1) * fwidth/2.;
627 else
628 posi -= foffset + (2*copyNo + 1) * fwidth/2.;
629
630 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
631 }
632
633 //----- calculate rotation matrix: unit
634
635#ifdef G4DIVDEBUG
636 if( verbose >= 2 )
637 {
638 G4cout << " G4ParameterisationPolyhedraZ - position: " << posi << G4endl
639 << " copyNo: " << copyNo << " - foffset: " << foffset/deg
640 << " - fwidth: " << fwidth/deg << G4endl;
641 }
642#endif
643
644 ChangeRotMatrix( physVol );
645
646#ifdef G4DIVDEBUG
647 if( verbose >= 2 )
648 {
649 G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraZ "
650 << copyNo << G4endl
651 << " Position: " << origin << " - Width: " << fwidth
652 << " - Axis: " << faxis << G4endl;
653 }
654#endif
655}
656
657//---------------------------------------------------------------------
658void
659G4ParameterisationPolyhedraZ::
660ComputeDimensions( G4Polyhedra& phedra, const G4int copyNo,
661 const G4VPhysicalVolume* ) const
662{
663 // Define division solid
664 G4PolyhedraHistorical origparam;
665 G4int nz = 2;
666 origparam.Num_z_planes = nz;
667 origparam.numSide = fOrigParamMother->numSide;
668 origparam.Start_angle = fOrigParamMother->Start_angle;
669 origparam.Opening_angle = fOrigParamMother->Opening_angle;
670
671 // Define division solid z sections
672 origparam.Z_values = new G4double[nz];
673 origparam.Rmin = new G4double[nz];
674 origparam.Rmax = new G4double[nz];
675 origparam.Z_values[0] = - fwidth/2.;
676 origparam.Z_values[1] = fwidth/2.;
677
678 if ( fDivisionType == DivNDIV ) {
679 // The position of the centre of copyNo-th mother polycone segment
680 G4double posi = ( fOrigParamMother->Z_values[copyNo]
681 + fOrigParamMother->Z_values[copyNo+1])/2;
682
683 origparam.Z_values[0] = fOrigParamMother->Z_values[copyNo] - posi;
684 origparam.Z_values[1] = fOrigParamMother->Z_values[copyNo+1] - posi;
685 origparam.Rmin[0] = fOrigParamMother->Rmin[copyNo];
686 origparam.Rmin[1] = fOrigParamMother->Rmin[copyNo+1];
687 origparam.Rmax[0] = fOrigParamMother->Rmax[copyNo];
688 origparam.Rmax[1] = fOrigParamMother->Rmax[copyNo+1];
689 }
690
691 if ( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH ) {
692 if ( ! fReflectedSolid ) {
693 origparam.Z_values[0] = - fwidth/2.;
694 origparam.Z_values[1] = fwidth/2.;
695
696 // The position of the centre of copyNo-th division
697 G4double posi
698 = fOrigParamMother->Z_values[0] + foffset + (2*copyNo + 1) * fwidth/2.;
699
700 // The first and last z sides z values
701 G4double zstart = posi - fwidth/2.;
702 G4double zend = posi + fwidth/2.;
703 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
704 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
705 origparam.Rmin[1] = GetRmin(zend, fNSegment);
706 origparam.Rmax[1] = GetRmax(zend, fNSegment);
707 }
708 else {
709 origparam.Z_values[0] = fwidth/2.;
710 origparam.Z_values[1] = - fwidth/2.;
711
712 // The position of the centre of copyNo-th division
713 G4double posi
714 = fOrigParamMother->Z_values[0] - ( foffset + (2*copyNo + 1) * fwidth/2.);
715
716 // The first and last z sides z values
717 G4double zstart = posi + fwidth/2.;
718 G4double zend = posi - fwidth/2.;
719 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
720 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
721 origparam.Rmin[1] = GetRmin(zend, fNSegment);
722 origparam.Rmax[1] = GetRmax(zend, fNSegment);
723 }
724
725 // It can happen due to rounding errors
726 if ( origparam.Rmin[0] < 0.0 ) origparam.Rmin[0] = 0.0;
727 if ( origparam.Rmin[nz-1] < 0.0 ) origparam.Rmin[1] = 0.0;
728 }
729
730 phedra.SetOriginalParameters(&origparam); // copy values & transfer pointers
731 phedra.Reset(); // reset to new solid parameters
732
733#ifdef G4DIVDEBUG
734 if( verbose >= 2 )
735 {
736 G4cout << "G4ParameterisationPolyhedraZ::ComputeDimensions()" << G4endl
737 << "-- Parametrised phedra copy-number: " << copyNo << G4endl;
738 phedra.DumpInfo();
739 }
740#endif
741}
742
Note: See TracBrowser for help on using the repository browser.