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