G4ParameterisationPolyhedra.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id: G4ParameterisationPolyhedra.cc 69784 2013-05-15 09:16:06Z gcosmo $
00028 //
00029 // class G4ParameterisationPolyhedra Implementation file
00030 //
00031 // 14.10.03 - P.Arce, Initial version
00032 // 08.04.04 - I.Hrivnacova, Implemented reflection
00033 // --------------------------------------------------------------------
00034 
00035 #include "G4ParameterisationPolyhedra.hh"
00036 
00037 #include <iomanip>
00038 #include "G4PhysicalConstants.hh"
00039 #include "G4ThreeVector.hh"
00040 #include "G4GeometryTolerance.hh"
00041 #include "G4RotationMatrix.hh"
00042 #include "G4VPhysicalVolume.hh"
00043 #include "G4LogicalVolume.hh"
00044 #include "G4ReflectedSolid.hh"
00045 #include "G4Polyhedra.hh"
00046 
00047 //--------------------------------------------------------------------------
00048 G4VParameterisationPolyhedra::
00049 G4VParameterisationPolyhedra( EAxis axis, G4int nDiv, G4double width,
00050                               G4double offset, G4VSolid* msolid,
00051                               DivisionType divType )
00052   :  G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid )
00053 {
00054   G4Polyhedra* msol = (G4Polyhedra*)(msolid);
00055   if ((msolid->GetEntityType() != "G4ReflectedSolid") && (msol->IsGeneric()))
00056   {
00057     std::ostringstream message;
00058     message << "Generic construct for G4Polyhedra NOT supported." << G4endl
00059             << "Sorry! Solid: " << msol->GetName();
00060     G4Exception("G4VParameterisationPolyhedra::G4VParameterisationPolyhedra()",
00061                 "GeomDiv0001", FatalException, message);
00062   }
00063   if (msolid->GetEntityType() == "G4ReflectedSolid")
00064   {
00065     // Get constituent solid  
00066     G4VSolid* mConstituentSolid 
00067        = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
00068     msol = (G4Polyhedra*)(mConstituentSolid);
00069   
00070     // Get parameters
00071     G4int   nofSides = msol->GetOriginalParameters()->numSide;
00072     G4int   nofZplanes = msol->GetOriginalParameters()->Num_z_planes;
00073     G4double* zValues  = msol->GetOriginalParameters()->Z_values;
00074     G4double* rminValues  = msol->GetOriginalParameters()->Rmin;
00075     G4double* rmaxValues  = msol->GetOriginalParameters()->Rmax;
00076 
00077     // Invert z values, 
00078     // convert radius parameters
00079     G4double* rminValues2 = new G4double[nofZplanes];
00080     G4double* rmaxValues2 = new G4double[nofZplanes];
00081     G4double* zValuesRefl = new G4double[nofZplanes];
00082     for (G4int i=0; i<nofZplanes; i++)
00083     {
00084       rminValues2[i] = rminValues[i] * ConvertRadiusFactor(*msol);
00085       rmaxValues2[i] = rmaxValues[i] * ConvertRadiusFactor(*msol);
00086       zValuesRefl[i] = - zValues[i];
00087     }  
00088     
00089     G4Polyhedra* newSolid
00090       = new G4Polyhedra(msol->GetName(),
00091                         msol->GetStartPhi(), 
00092                         msol->GetEndPhi() - msol->GetStartPhi(),
00093                         nofSides,
00094                         nofZplanes, zValuesRefl, rminValues2, rmaxValues2);
00095 
00096     delete [] rminValues2;       
00097     delete [] rmaxValues2;       
00098     delete [] zValuesRefl;       
00099 
00100     msol = newSolid;
00101     fmotherSolid = newSolid;
00102     fReflectedSolid = true;
00103     fDeleteSolid = true;
00104   }    
00105 }
00106 
00107 //------------------------------------------------------------------------
00108 G4VParameterisationPolyhedra::~G4VParameterisationPolyhedra()
00109 {
00110 }
00111 
00112 //--------------------------------------------------------------------------
00113 G4double 
00114 G4VParameterisationPolyhedra::
00115 ConvertRadiusFactor(const G4Polyhedra& phedra) const
00116 {
00117   G4double phiTotal = phedra.GetEndPhi() - phedra.GetStartPhi();
00118   G4int nofSides = phedra.GetOriginalParameters()->numSide;
00119   
00120   if ( (phiTotal <=0) || (phiTotal >
00121         2*pi+G4GeometryTolerance::GetInstance()->GetAngularTolerance()) )
00122     { phiTotal = 2*pi; }
00123   
00124   return std::cos(0.5*phiTotal/nofSides);
00125 }  
00126 
00127 //--------------------------------------------------------------------------
00128 G4ParameterisationPolyhedraRho::
00129 G4ParameterisationPolyhedraRho( EAxis axis, G4int nDiv,
00130                                G4double width, G4double offset,
00131                                G4VSolid* msolid, DivisionType divType )
00132   :  G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType )
00133 {
00134 
00135   CheckParametersValidity();
00136   SetType( "DivisionPolyhedraRho" );
00137 
00138   G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
00139   G4PolyhedraHistorical* original_pars = msol->GetOriginalParameters();
00140 
00141   if( divType == DivWIDTH )
00142   {
00143     fnDiv = CalculateNDiv( original_pars->Rmax[0]
00144                          - original_pars->Rmin[0], width, offset );
00145   }
00146   else if( divType == DivNDIV )
00147   {
00148     fwidth = CalculateWidth( original_pars->Rmax[0]
00149                            - original_pars->Rmin[0], nDiv, offset );
00150   }
00151 
00152 #ifdef G4DIVDEBUG
00153   if( verbose >= 1 )
00154   {
00155     G4cout << " G4ParameterisationPolyhedraRho - # divisions " << fnDiv
00156            << " = " << nDiv << G4endl
00157            << " Offset " << foffset << " = " << offset << G4endl
00158            << " Width " << fwidth << " = " << width << G4endl;
00159   }
00160 #endif
00161 }
00162 
00163 //------------------------------------------------------------------------
00164 G4ParameterisationPolyhedraRho::~G4ParameterisationPolyhedraRho()
00165 {
00166 }
00167 
00168 //---------------------------------------------------------------------
00169 void G4ParameterisationPolyhedraRho::CheckParametersValidity()
00170 {
00171   G4VDivisionParameterisation::CheckParametersValidity();
00172 
00173   G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
00174 
00175   if( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH )
00176   {
00177     std::ostringstream message;
00178     message << "In solid " << msol->GetName() << G4endl
00179             << "Division along R will be done with a width "
00180             << "different for each solid section." << G4endl
00181             << "WIDTH will not be used !";
00182     G4Exception("G4ParameterisationPolyhedraRho::CheckParametersValidity()",
00183                 "GeomDiv1001", JustWarning, message);
00184   }
00185   if( foffset != 0. )
00186   {
00187     std::ostringstream message;
00188     message << "In solid " << msol->GetName() << G4endl
00189             << "Division along  R will be done with a width "
00190             << "different for each solid section." << G4endl
00191             << "OFFSET will not be used !";
00192     G4Exception("G4ParameterisationPolyhedraRho::CheckParametersValidity()",
00193                 "GeomDiv1001", JustWarning, message);
00194   }
00195 }
00196 
00197 //------------------------------------------------------------------------
00198 G4double G4ParameterisationPolyhedraRho::GetMaxParameter() const
00199 {
00200   G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
00201   G4PolyhedraHistorical* original_pars = msol->GetOriginalParameters();
00202   return original_pars->Rmax[0] - original_pars->Rmin[0];
00203 }
00204 
00205 //--------------------------------------------------------------------------
00206 void
00207 G4ParameterisationPolyhedraRho::
00208 ComputeTransformation( const G4int, G4VPhysicalVolume* physVol ) const
00209 {
00210   //----- translation 
00211   G4ThreeVector origin(0.,0.,0.);
00212 
00213   //----- set translation 
00214   physVol->SetTranslation( origin );
00215 
00216   //----- calculate rotation matrix: unit
00217 
00218 #ifdef G4DIVDEBUG
00219   if( verbose >= 2 )
00220   {
00221     G4cout << " G4ParameterisationPolyhedraRho " << G4endl
00222            << " foffset: " << foffset/deg
00223            << " - fwidth: " << fwidth/deg << G4endl;
00224   }
00225 #endif
00226 
00227   ChangeRotMatrix( physVol );
00228 
00229 #ifdef G4DIVDEBUG
00230   if( verbose >= 2 )
00231   {
00232     G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraRho "
00233            << G4endl
00234            << " Position: " << origin
00235            << " - Width: " << fwidth
00236            << " - Axis: " << faxis  << G4endl;
00237   }
00238 #endif
00239 }
00240 
00241 //--------------------------------------------------------------------------
00242 void
00243 G4ParameterisationPolyhedraRho::
00244 ComputeDimensions( G4Polyhedra& phedra, const G4int copyNo,
00245                    const G4VPhysicalVolume* ) const
00246 {
00247   G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
00248 
00249   G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters();
00250   G4PolyhedraHistorical origparam( *origparamMother );
00251   G4int nZplanes = origparamMother->Num_z_planes;
00252 
00253   G4double width = 0.;
00254   for( G4int ii = 0; ii < nZplanes; ii++ )
00255   {
00256     width = CalculateWidth( origparamMother->Rmax[ii]
00257                           - origparamMother->Rmin[ii], fnDiv, foffset );
00258     origparam.Rmin[ii] = origparamMother->Rmin[ii]+foffset+width*copyNo;
00259     origparam.Rmax[ii] = origparamMother->Rmin[ii]+foffset+width*(copyNo+1);
00260   }
00261 
00262   phedra.SetOriginalParameters(&origparam); // copy values & transfer pointers
00263   phedra.Reset();                           // reset to new solid parameters
00264 
00265 #ifdef G4DIVDEBUG
00266   if( verbose >= -2 )
00267   {
00268     G4cout << "G4ParameterisationPolyhedraRho::ComputeDimensions()" << G4endl
00269            << "-- Parametrised phedra copy-number: " << copyNo << G4endl;
00270     phedra.DumpInfo();
00271   }
00272 #endif
00273 }
00274 
00275 //--------------------------------------------------------------------------
00276 G4ParameterisationPolyhedraPhi::
00277 G4ParameterisationPolyhedraPhi( EAxis axis, G4int nDiv,
00278                                G4double width, G4double offset,
00279                                G4VSolid* msolid, DivisionType divType )
00280   :  G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType )
00281 { 
00282   CheckParametersValidity();
00283   SetType( "DivisionPolyhedraPhi" );
00284 
00285   G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
00286   G4double deltaPhi = msol->GetEndPhi() - msol->GetStartPhi();
00287 
00288   if( divType == DivWIDTH )
00289   {
00290     fnDiv = msol->GetNumSide();
00291   }
00292 
00293   fwidth = CalculateWidth( deltaPhi, fnDiv, 0.0 );
00294 
00295 #ifdef G4DIVDEBUG
00296   if( verbose >= 1 )
00297   {
00298     G4cout << " G4ParameterisationPolyhedraPhi - # divisions " << fnDiv
00299            << " = " << nDiv << G4endl
00300            << " Offset " << foffset << " = " << offset << G4endl
00301            << " Width " << fwidth << " = " << width << G4endl;
00302   }
00303 #endif
00304 }
00305 
00306 //------------------------------------------------------------------------
00307 G4ParameterisationPolyhedraPhi::~G4ParameterisationPolyhedraPhi()
00308 {
00309 }
00310 
00311 //------------------------------------------------------------------------
00312 G4double G4ParameterisationPolyhedraPhi::GetMaxParameter() const
00313 {
00314   G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
00315   return msol->GetEndPhi() - msol->GetStartPhi();
00316 }
00317 
00318 //---------------------------------------------------------------------
00319 void G4ParameterisationPolyhedraPhi::CheckParametersValidity()
00320 {
00321   G4VDivisionParameterisation::CheckParametersValidity();
00322 
00323   G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
00324 
00325   if( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH )
00326   {
00327     std::ostringstream message;
00328     message << "In solid " << msol->GetName() << G4endl
00329             << " Division along PHI will be done splitting "
00330             << "in the defined numSide." << G4endl
00331             << "WIDTH will not be used !";
00332     G4Exception("G4ParameterisationPolyhedraPhi::CheckParametersValidity()",
00333                 "GeomDiv1001", JustWarning, message);
00334   }
00335   if( foffset != 0. )
00336   {
00337     std::ostringstream message;
00338     message << "In solid " << msol->GetName() << G4endl
00339             << "Division along PHI will be done splitting "
00340             << "in the defined numSide." << G4endl
00341             << "OFFSET will not be used !";
00342     G4Exception("G4ParameterisationPolyhedraPhi::CheckParametersValidity()",
00343                 "GeomDiv1001", JustWarning, message);
00344   }
00345 
00346   G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters();
00347 
00348   if( origparamMother->numSide != fnDiv &&  fDivisionType != DivWIDTH)
00349   { 
00350     std::ostringstream message;
00351     message << "Configuration not supported." << G4endl
00352             << "Division along PHI will be done splitting in the defined"
00353             << G4endl
00354             << "numSide, i.e, the number of division would be :"
00355             << origparamMother->numSide << " instead of " << fnDiv << " !"; 
00356     G4Exception("G4ParameterisationPolyhedraPhi::CheckParametersValidity()",
00357                 "GeomDiv0001", FatalException, message);
00358   }
00359 }
00360 
00361 //--------------------------------------------------------------------------
00362 void
00363 G4ParameterisationPolyhedraPhi::
00364 ComputeTransformation( const G4int copyNo, G4VPhysicalVolume *physVol ) const
00365 {
00366   //----- translation 
00367   G4ThreeVector origin(0.,0.,0.); 
00368   //----- set translation 
00369   physVol->SetTranslation( origin );
00370 
00371   //----- calculate rotation matrix (so that all volumes point to the centre)
00372   G4double posi = copyNo*fwidth;
00373 
00374 #ifdef G4DIVDEBUG
00375   if( verbose >= 2 )
00376   {
00377     G4cout << " G4ParameterisationPolyhedraPhi - position: " << posi/deg
00378            << G4endl
00379            << " copyNo: " << copyNo
00380            << " - fwidth: " << fwidth/deg << G4endl;
00381   }
00382 #endif
00383 
00384   ChangeRotMatrix( physVol, -posi );
00385 
00386 #ifdef G4DIVDEBUG
00387   if( verbose >= 2 )
00388   {
00389     G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraPhi " << copyNo
00390            << G4endl
00391            << " Position: " << origin << " - Width: " << fwidth
00392            << " - Axis: " << faxis  << G4endl;
00393   }
00394 #endif
00395 }
00396 
00397 //--------------------------------------------------------------------------
00398 void
00399 G4ParameterisationPolyhedraPhi::
00400 ComputeDimensions( G4Polyhedra& phedra, const G4int,
00401                    const G4VPhysicalVolume* ) const
00402 {
00403   G4Polyhedra* msol = (G4Polyhedra*)(fmotherSolid);
00404 
00405   G4PolyhedraHistorical* origparamMother = msol->GetOriginalParameters();
00406   G4PolyhedraHistorical origparam( *origparamMother );
00407 
00408   origparam.numSide = 1;
00409   origparam.Start_angle = origparamMother->Start_angle;
00410   origparam.Opening_angle = fwidth;
00411 
00412   phedra.SetOriginalParameters(&origparam);  // copy values & transfer pointers
00413   phedra.Reset();                            // reset to new solid parameters
00414 
00415 #ifdef G4DIVDEBUG
00416   if( verbose >= 2 )
00417   {
00418     G4cout << "G4ParameterisationPolyhedraPhi::ComputeDimensions():" << G4endl;
00419     phedra.DumpInfo();
00420   }
00421 #endif
00422 }
00423 
00424 //--------------------------------------------------------------------------
00425 G4ParameterisationPolyhedraZ::
00426 G4ParameterisationPolyhedraZ( EAxis axis, G4int nDiv,
00427                              G4double width, G4double offset,
00428                              G4VSolid* msolid, DivisionType divType )
00429   :  G4VParameterisationPolyhedra( axis, nDiv, width, offset, msolid, divType ),
00430      fNSegment(0),
00431      fOrigParamMother(((G4Polyhedra*)fmotherSolid)->GetOriginalParameters())
00432 { 
00433   CheckParametersValidity();
00434   SetType( "DivisionPolyhedraZ" );
00435 
00436   if( divType == DivWIDTH )
00437     {
00438     fnDiv =
00439       CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
00440                      - fOrigParamMother->Z_values[0] , width, offset );
00441   }
00442   else if( divType == DivNDIV )
00443     {
00444     fwidth =
00445       CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
00446                      - fOrigParamMother->Z_values[0] , nDiv, offset );
00447   }
00448   
00449 #ifdef G4DIVDEBUG
00450   if( verbose >= 1 )
00451     {
00452     G4cout << " G4ParameterisationPolyhedraZ - # divisions " << fnDiv << " = "
00453            << nDiv << G4endl
00454            << " Offset " << foffset << " = " << offset << G4endl
00455            << " Width " << fwidth << " = " << width << G4endl;
00456   }
00457 #endif
00458 }
00459 
00460 //---------------------------------------------------------------------
00461 G4ParameterisationPolyhedraZ::~G4ParameterisationPolyhedraZ()
00462 {
00463 }
00464 
00465 //------------------------------------------------------------------------
00466 G4double G4ParameterisationPolyhedraZ::GetR(G4double z, 
00467                                            G4double z1, G4double r1,  
00468                                            G4double z2, G4double r2) const
00469 {
00470   // Linear parameterisation: 
00471   // r = az + b
00472   // a = (r1 - r2)/(z1-z2)
00473   // b = r1 - a*z1
00474 
00475   return (r1-r2)/(z1-z2)*z + ( r1 - (r1-r2)/(z1-z2)*z1 ) ;
00476 }  
00477                                            
00478 //------------------------------------------------------------------------
00479 G4double G4ParameterisationPolyhedraZ::GetRmin(G4double z, G4int nseg) const
00480 {
00481 // Get Rmin in the given z position for the given polyhedra segment 
00482 
00483   return GetR(z, 
00484               fOrigParamMother->Z_values[nseg], 
00485               fOrigParamMother->Rmin[nseg],
00486               fOrigParamMother->Z_values[nseg+1], 
00487               fOrigParamMother->Rmin[nseg+1]);
00488 }  
00489                                            
00490 //------------------------------------------------------------------------
00491 G4double G4ParameterisationPolyhedraZ::GetRmax(G4double z, G4int nseg) const
00492 {
00493 // Get Rmax in the given z position for the given polyhedra segment 
00494 
00495   return GetR(z, 
00496               fOrigParamMother->Z_values[nseg], 
00497               fOrigParamMother->Rmax[nseg],
00498               fOrigParamMother->Z_values[nseg+1], 
00499               fOrigParamMother->Rmax[nseg+1]);
00500 }  
00501                                            
00502 //------------------------------------------------------------------------
00503 G4double G4ParameterisationPolyhedraZ::GetMaxParameter() const
00504 {
00505   return std::abs (fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
00506              -fOrigParamMother->Z_values[0]);
00507 }
00508 
00509 //---------------------------------------------------------------------
00510 void G4ParameterisationPolyhedraZ::CheckParametersValidity()
00511 {
00512   G4VDivisionParameterisation::CheckParametersValidity();
00513 
00514   // Division will be following the mother polyhedra segments
00515   if( fDivisionType == DivNDIV ) {
00516     if( fOrigParamMother->Num_z_planes-1 != fnDiv ) { 
00517       std::ostringstream message;
00518       message << "Configuration not supported." << G4endl
00519               << "Division along Z will be done splitting in the defined"
00520               << G4endl
00521               << "Z planes, i.e, the number of division would be :"
00522               << fOrigParamMother->Num_z_planes-1 << " instead of "
00523               << fnDiv << " !"; 
00524       G4Exception("G4ParameterisationPolyhedraZ::CheckParametersValidity()",
00525                   "GeomDiv0001", FatalException, message);
00526     }
00527   }  
00528 
00529   // Division will be done within one polyhedra segment
00530   // with applying given width and offset
00531   if( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH ) {
00532     // Check if divided region does not span over more
00533     // than one z segment
00534 
00535     G4int isegstart = -1;  // number of the segment containing start position
00536     G4int isegend = -1;    // number of the segment containing end position
00537 
00538     if ( ! fReflectedSolid ) {
00539       // The start/end position of the divided region
00540       G4double zstart 
00541         = fOrigParamMother->Z_values[0] + foffset;
00542       G4double zend 
00543         = fOrigParamMother->Z_values[0] + foffset + fnDiv* fwidth;
00544    
00545       G4int counter = 0;
00546       while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 ) {
00547         // first segment
00548         if ( zstart >= fOrigParamMother->Z_values[counter]  &&
00549              zstart  < fOrigParamMother->Z_values[counter+1] ) {
00550            isegstart = counter;
00551         }     
00552         // last segment
00553         if ( zend  > fOrigParamMother->Z_values[counter] &&
00554              zend <= fOrigParamMother->Z_values[counter+1] ) {
00555            isegend = counter;
00556         }   
00557         ++counter;   
00558       }
00559     }
00560     else  {
00561       // The start/end position of the divided region
00562       G4double zstart 
00563         = fOrigParamMother->Z_values[0] - foffset;
00564       G4double zend 
00565         = fOrigParamMother->Z_values[0] - ( foffset + fnDiv* fwidth);
00566    
00567       G4int counter = 0;
00568       while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 ) {
00569         // first segment
00570         if ( zstart <= fOrigParamMother->Z_values[counter]  &&
00571              zstart  > fOrigParamMother->Z_values[counter+1] ) {
00572            isegstart = counter;
00573         }     
00574         // last segment
00575         if ( zend  < fOrigParamMother->Z_values[counter] &&
00576              zend >= fOrigParamMother->Z_values[counter+1] ) {
00577            isegend = counter;
00578         }   
00579         ++counter;   
00580       }
00581     }
00582   
00583     if ( isegstart != isegend ) {
00584       std::ostringstream message;
00585       message << "Configuration not supported." << G4endl
00586               << "Division with user defined width." << G4endl
00587               << "Solid " << fmotherSolid->GetName() << G4endl
00588               << "Divided region is not between two Z planes."; 
00589       G4Exception("G4ParameterisationPolyhedraZ::CheckParametersValidity()",
00590                   "GeomDiv0001", FatalException, message);
00591     }
00592   
00593     fNSegment = isegstart;
00594   }  
00595 }
00596 
00597 //---------------------------------------------------------------------
00598 void
00599 G4ParameterisationPolyhedraZ::
00600 ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol) const
00601 {
00602   if ( fDivisionType == DivNDIV ) {
00603     // The position of the centre of copyNo-th mother polycone segment
00604     G4double posi = ( fOrigParamMother->Z_values[copyNo]
00605                     + fOrigParamMother->Z_values[copyNo+1])/2;
00606     physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
00607   }
00608   
00609   if ( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH ) {
00610     // The position of the centre of copyNo-th division
00611 
00612     G4double posi = fOrigParamMother->Z_values[0];
00613     
00614     if ( ! fReflectedSolid )
00615       posi += foffset + (2*copyNo + 1) * fwidth/2.;
00616     else
00617       posi -= foffset + (2*copyNo + 1) * fwidth/2.;
00618     
00619     physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
00620   }   
00621 
00622   //----- calculate rotation matrix: unit
00623 
00624 #ifdef G4DIVDEBUG
00625   if( verbose >= 2 )
00626   {
00627     G4cout << " G4ParameterisationPolyhedraZ - position: " << posi << G4endl
00628            << " copyNo: " << copyNo << " - foffset: " << foffset/deg
00629            << " - fwidth: " << fwidth/deg << G4endl;
00630   }
00631 #endif
00632 
00633   ChangeRotMatrix( physVol );
00634 
00635 #ifdef G4DIVDEBUG
00636   if( verbose >= 2 )
00637   {
00638     G4cout << std::setprecision(8) << " G4ParameterisationPolyhedraZ "
00639            << copyNo << G4endl
00640            << " Position: " << origin << " - Width: " << fwidth
00641            << " - Axis: " << faxis  << G4endl;
00642   }
00643 #endif
00644 }
00645 
00646 //---------------------------------------------------------------------
00647 void
00648 G4ParameterisationPolyhedraZ::
00649 ComputeDimensions( G4Polyhedra& phedra, const G4int copyNo,
00650                    const G4VPhysicalVolume* ) const
00651 {
00652   // Define division solid
00653   G4PolyhedraHistorical origparam;
00654   G4int nz = 2; 
00655   origparam.Num_z_planes = nz;
00656   origparam.numSide = fOrigParamMother->numSide;
00657   origparam.Start_angle = fOrigParamMother->Start_angle;
00658   origparam.Opening_angle = fOrigParamMother->Opening_angle;
00659 
00660   // Define division solid z sections
00661   origparam.Z_values = new G4double[nz];
00662   origparam.Rmin = new G4double[nz];
00663   origparam.Rmax = new G4double[nz];
00664   origparam.Z_values[0] = - fwidth/2.;
00665   origparam.Z_values[1] = fwidth/2.;
00666 
00667   if ( fDivisionType == DivNDIV ) {
00668     // The position of the centre of copyNo-th mother polycone segment
00669     G4double posi = ( fOrigParamMother->Z_values[copyNo]
00670                     + fOrigParamMother->Z_values[copyNo+1])/2;
00671 
00672     origparam.Z_values[0] = fOrigParamMother->Z_values[copyNo] - posi;
00673     origparam.Z_values[1] = fOrigParamMother->Z_values[copyNo+1] - posi;
00674     origparam.Rmin[0] = fOrigParamMother->Rmin[copyNo];
00675     origparam.Rmin[1] = fOrigParamMother->Rmin[copyNo+1];
00676     origparam.Rmax[0] = fOrigParamMother->Rmax[copyNo];
00677     origparam.Rmax[1] = fOrigParamMother->Rmax[copyNo+1];
00678   }  
00679 
00680   if ( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH ) {
00681     if ( ! fReflectedSolid ) {
00682       origparam.Z_values[0] = - fwidth/2.;
00683       origparam.Z_values[1] = fwidth/2.;
00684 
00685       // The position of the centre of copyNo-th division
00686       G4double posi 
00687         = fOrigParamMother->Z_values[0] + foffset + (2*copyNo + 1) * fwidth/2.;
00688     
00689       // The first and last z sides z values
00690       G4double zstart = posi - fwidth/2.;
00691       G4double zend = posi + fwidth/2.;
00692       origparam.Rmin[0] = GetRmin(zstart, fNSegment); 
00693       origparam.Rmax[0] = GetRmax(zstart, fNSegment);  
00694       origparam.Rmin[1] = GetRmin(zend, fNSegment); 
00695       origparam.Rmax[1] = GetRmax(zend, fNSegment); 
00696     }
00697     else {   
00698       origparam.Z_values[0] = fwidth/2.;
00699       origparam.Z_values[1] = - fwidth/2.;
00700 
00701       // The position of the centre of copyNo-th division
00702       G4double posi 
00703         = fOrigParamMother->Z_values[0] - ( foffset + (2*copyNo + 1) * fwidth/2.);
00704     
00705       // The first and last z sides z values
00706       G4double zstart = posi + fwidth/2.;
00707       G4double zend = posi - fwidth/2.;
00708       origparam.Rmin[0] = GetRmin(zstart, fNSegment); 
00709       origparam.Rmax[0] = GetRmax(zstart, fNSegment);  
00710       origparam.Rmin[1] = GetRmin(zend, fNSegment); 
00711       origparam.Rmax[1] = GetRmax(zend, fNSegment); 
00712     } 
00713 
00714     // It can happen due to rounding errors
00715     if ( origparam.Rmin[0]    < 0.0 ) origparam.Rmin[0] = 0.0;
00716     if ( origparam.Rmin[nz-1] < 0.0 ) origparam.Rmin[1] = 0.0;
00717   }  
00718 
00719   phedra.SetOriginalParameters(&origparam);  // copy values & transfer pointers
00720   phedra.Reset();                            // reset to new solid parameters
00721 
00722 #ifdef G4DIVDEBUG
00723   if( verbose >= 2 )
00724   {
00725     G4cout << "G4ParameterisationPolyhedraZ::ComputeDimensions()" << G4endl
00726            << "-- Parametrised phedra copy-number: " << copyNo << G4endl;
00727     phedra.DumpInfo();
00728   }
00729 #endif
00730 }

Generated on Mon May 27 17:49:15 2013 for Geant4 by  doxygen 1.4.7