00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035 #include "G4ParameterisationPolycone.hh"
00036
00037 #include <iomanip>
00038 #include "G4ThreeVector.hh"
00039 #include "G4RotationMatrix.hh"
00040 #include "G4VPhysicalVolume.hh"
00041 #include "G4LogicalVolume.hh"
00042 #include "G4ReflectedSolid.hh"
00043
00044
00045 G4VParameterisationPolycone::
00046 G4VParameterisationPolycone( EAxis axis, G4int nDiv, G4double width,
00047 G4double offset, G4VSolid* msolid,
00048 DivisionType divType )
00049 : G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid )
00050 {
00051 G4Polycone* msol = (G4Polycone*)(msolid);
00052 if ((msolid->GetEntityType() != "G4ReflectedSolid") && (msol->IsGeneric()))
00053 {
00054 std::ostringstream message;
00055 message << "Generic construct for G4Polycone NOT supported." << G4endl
00056 << "Sorry! Solid: " << msol->GetName() << ".";
00057 G4Exception("G4VParameterisationPolycone::G4VParameterisationPolycone()",
00058 "GeomDiv0001", FatalException, message);
00059 }
00060 if (msolid->GetEntityType() == "G4ReflectedSolid")
00061 {
00062
00063 G4VSolid* mConstituentSolid
00064 = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
00065 msol = (G4Polycone*)(mConstituentSolid);
00066
00067
00068 G4int nofZplanes = msol->GetOriginalParameters()->Num_z_planes;
00069 G4double* zValues = msol->GetOriginalParameters()->Z_values;
00070 G4double* rminValues = msol->GetOriginalParameters()->Rmin;
00071 G4double* rmaxValues = msol->GetOriginalParameters()->Rmax;
00072
00073
00074 G4double* zValuesRefl = new double[nofZplanes];
00075 for (G4int i=0; i<nofZplanes; i++) zValuesRefl[i] = - zValues[i];
00076
00077 G4Polycone* newSolid
00078 = new G4Polycone(msol->GetName(),
00079 msol->GetStartPhi(),
00080 msol->GetEndPhi() - msol->GetStartPhi(),
00081 nofZplanes, zValuesRefl, rminValues, rmaxValues);
00082
00083 delete [] zValuesRefl;
00084
00085 msol = newSolid;
00086 fmotherSolid = newSolid;
00087 fReflectedSolid = true;
00088 fDeleteSolid = true;
00089 }
00090 }
00091
00092
00093 G4VParameterisationPolycone::~G4VParameterisationPolycone()
00094 {
00095 }
00096
00097
00098 G4ParameterisationPolyconeRho::
00099 G4ParameterisationPolyconeRho( EAxis axis, G4int nDiv,
00100 G4double width, G4double offset,
00101 G4VSolid* msolid, DivisionType divType )
00102 : G4VParameterisationPolycone( axis, nDiv, width, offset, msolid, divType )
00103 {
00104 CheckParametersValidity();
00105 SetType( "DivisionPolyconeRho" );
00106
00107 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
00108 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
00109
00110 if( divType == DivWIDTH )
00111 {
00112 fnDiv = CalculateNDiv( origparamMother->Rmax[0]
00113 - origparamMother->Rmin[0], width, offset );
00114 }
00115 else if( divType == DivNDIV )
00116 {
00117 fwidth = CalculateWidth( origparamMother->Rmax[0]
00118 - origparamMother->Rmin[0], nDiv, offset );
00119 }
00120
00121 #ifdef G4DIVDEBUG
00122 if( verbose >= -1 )
00123 {
00124 G4cout << " G4ParameterisationPolyconeRho - # divisions " << fnDiv
00125 << " = " << nDiv << G4endl
00126 << " Offset " << foffset << " = " << offset << G4endl
00127 << " Width " << fwidth << " = " << width << G4endl;
00128 }
00129 #endif
00130 }
00131
00132
00133 G4ParameterisationPolyconeRho::~G4ParameterisationPolyconeRho()
00134 {
00135 }
00136
00137
00138 void G4ParameterisationPolyconeRho::CheckParametersValidity()
00139 {
00140 G4VDivisionParameterisation::CheckParametersValidity();
00141
00142 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
00143
00144 if( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH )
00145 {
00146 std::ostringstream message;
00147 message << "In solid " << msol->GetName() << G4endl
00148 << "Division along R will be done with a width "
00149 << "different for each solid section." << G4endl
00150 << "WIDTH will not be used !";
00151 G4Exception("G4VParameterisationPolycone::CheckParametersValidity()",
00152 "GeomDiv1001", JustWarning, message);
00153 }
00154 if( foffset != 0. )
00155 {
00156 std::ostringstream message;
00157 message << "In solid " << msol->GetName() << G4endl
00158 << "Division along R will be done with a width "
00159 << "different for each solid section." << G4endl
00160 << "OFFSET will not be used !";
00161 G4Exception("G4VParameterisationPolycone::CheckParametersValidity()",
00162 "GeomDiv1001", JustWarning, message);
00163 }
00164 }
00165
00166
00167 G4double G4ParameterisationPolyconeRho::GetMaxParameter() const
00168 {
00169 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
00170 G4PolyconeHistorical* original_pars = msol->GetOriginalParameters();
00171 return original_pars->Rmax[0] - original_pars->Rmin[0];
00172 }
00173
00174
00175
00176 void
00177 G4ParameterisationPolyconeRho::
00178 ComputeTransformation( const G4int, G4VPhysicalVolume* physVol ) const
00179 {
00180
00181 G4ThreeVector origin(0.,0.,0.);
00182
00183 physVol->SetTranslation( origin );
00184
00185
00186
00187 #ifdef G4DIVDEBUG
00188 if( verbose >= 2 )
00189 {
00190 G4cout << " G4ParameterisationPolyconeRho " << G4endl
00191 << " foffset: " << foffset
00192 << " - fwidth: " << fwidth << G4endl;
00193 }
00194 #endif
00195
00196 ChangeRotMatrix( physVol );
00197
00198 #ifdef G4DIVDEBUG
00199 if( verbose >= 2 )
00200 {
00201 G4cout << std::setprecision(8) << " G4ParameterisationPolyconeRho "
00202 << G4endl
00203 << " Position: " << origin/mm
00204 << " - Width: " << fwidth/deg
00205 << " - Axis: " << faxis << G4endl;
00206 }
00207 #endif
00208 }
00209
00210
00211 void
00212 G4ParameterisationPolyconeRho::
00213 ComputeDimensions( G4Polycone& pcone, const G4int copyNo,
00214 const G4VPhysicalVolume* ) const
00215 {
00216 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
00217
00218 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
00219 G4PolyconeHistorical origparam( *origparamMother );
00220 G4int nZplanes = origparamMother->Num_z_planes;
00221
00222 G4double width = 0.;
00223 for( G4int ii = 0; ii < nZplanes; ii++ )
00224 {
00225 width = CalculateWidth( origparamMother->Rmax[ii]
00226 - origparamMother->Rmin[ii], fnDiv, foffset );
00227 origparam.Rmin[ii] = origparamMother->Rmin[ii]+foffset+width*copyNo;
00228 origparam.Rmax[ii] = origparamMother->Rmin[ii]+foffset+width*(copyNo+1);
00229 }
00230
00231 pcone.SetOriginalParameters(&origparam);
00232 pcone.Reset();
00233
00234 #ifdef G4DIVDEBUG
00235 if( verbose >= -2 )
00236 {
00237 G4cout << "G4ParameterisationPolyconeRho::ComputeDimensions()" << G4endl
00238 << "-- Parametrised pcone copy-number: " << copyNo << G4endl;
00239 pcone.DumpInfo();
00240 }
00241 #endif
00242 }
00243
00244
00245 G4ParameterisationPolyconePhi::
00246 G4ParameterisationPolyconePhi( EAxis axis, G4int nDiv,
00247 G4double width, G4double offset,
00248 G4VSolid* msolid, DivisionType divType )
00249 : G4VParameterisationPolycone( axis, nDiv, width, offset, msolid, divType )
00250 {
00251 CheckParametersValidity();
00252 SetType( "DivisionPolyconePhi" );
00253
00254 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
00255 G4double deltaPhi = msol->GetEndPhi() - msol->GetStartPhi();
00256
00257 if( divType == DivWIDTH )
00258 {
00259 fnDiv = CalculateNDiv( deltaPhi, width, offset );
00260 }
00261 else if( divType == DivNDIV )
00262 {
00263 fwidth = CalculateWidth( deltaPhi, nDiv, offset );
00264 }
00265
00266 #ifdef G4DIVDEBUG
00267 if( verbose >= 1 )
00268 {
00269 G4cout << " G4ParameterisationPolyconePhi - # divisions " << fnDiv
00270 << " = " << nDiv << G4endl
00271 << " Offset " << foffset/deg << " = " << offset/deg << G4endl
00272 << " Width " << fwidth/deg << " = " << width/deg << G4endl;
00273 }
00274 #endif
00275 }
00276
00277
00278 G4ParameterisationPolyconePhi::~G4ParameterisationPolyconePhi()
00279 {
00280 }
00281
00282
00283 G4double G4ParameterisationPolyconePhi::GetMaxParameter() const
00284 {
00285 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
00286 return msol->GetEndPhi() - msol->GetStartPhi();
00287 }
00288
00289
00290 void
00291 G4ParameterisationPolyconePhi::
00292 ComputeTransformation( const G4int copyNo, G4VPhysicalVolume *physVol ) const
00293 {
00294
00295 G4ThreeVector origin(0.,0.,0.);
00296
00297 physVol->SetTranslation( origin );
00298
00299
00300 G4double posi = foffset + copyNo*fwidth;
00301
00302 #ifdef G4DIVDEBUG
00303 if( verbose >= 2 )
00304 {
00305 G4cout << " G4ParameterisationPolyconePhi - position: " << posi/deg
00306 << G4endl
00307 << " copyNo: " << copyNo << " - foffset: " << foffset/deg
00308 << " - fwidth: " << fwidth/deg << G4endl;
00309 }
00310 #endif
00311
00312 ChangeRotMatrix( physVol, -posi );
00313
00314 #ifdef G4DIVDEBUG
00315 if( verbose >= 2 )
00316 {
00317 G4cout << std::setprecision(8) << " G4ParameterisationPolyconePhi "
00318 << copyNo << G4endl
00319 << " Position: " << origin << " - Width: " << fwidth
00320 << " - Axis: " << faxis << G4endl;
00321 }
00322 #endif
00323 }
00324
00325
00326 void
00327 G4ParameterisationPolyconePhi::
00328 ComputeDimensions( G4Polycone& pcone, const G4int,
00329 const G4VPhysicalVolume* ) const
00330 {
00331 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
00332
00333 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
00334 G4PolyconeHistorical origparam( *origparamMother );
00335 origparam.Start_angle = origparamMother->Start_angle;
00336 origparam.Opening_angle = fwidth;
00337
00338 pcone.SetOriginalParameters(&origparam);
00339 pcone.Reset();
00340
00341 #ifdef G4DIVDEBUG
00342 if( verbose >= 2 )
00343 {
00344 G4cout << "G4ParameterisationPolyconePhi::ComputeDimensions():" << G4endl;
00345 pcone.DumpInfo();
00346 }
00347 #endif
00348 }
00349
00350
00351 G4ParameterisationPolyconeZ::
00352 G4ParameterisationPolyconeZ( EAxis axis, G4int nDiv,
00353 G4double width, G4double offset,
00354 G4VSolid* msolid, DivisionType divType)
00355 : G4VParameterisationPolycone( axis, nDiv, width, offset, msolid, divType ),
00356 fNSegment(0),
00357 fOrigParamMother(((G4Polycone*)fmotherSolid)->GetOriginalParameters())
00358 {
00359
00360 CheckParametersValidity();
00361 SetType( "DivisionPolyconeZ" );
00362
00363 if( divType == DivWIDTH )
00364 {
00365 fnDiv =
00366 CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
00367 - fOrigParamMother->Z_values[0] , width, offset );
00368 }
00369 else if( divType == DivNDIV )
00370 {
00371 fwidth =
00372 CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
00373 - fOrigParamMother->Z_values[0] , nDiv, offset );
00374 }
00375
00376 #ifdef G4DIVDEBUG
00377 if( verbose >= 1 )
00378 {
00379 G4cout << " G4ParameterisationPolyconeZ - # divisions " << fnDiv << " = "
00380 << nDiv << G4endl
00381 << " Offset " << foffset << " = " << offset << G4endl
00382 << " Width " << fwidth << " = " << width << G4endl;
00383 }
00384 #endif
00385 }
00386
00387
00388 G4ParameterisationPolyconeZ::~G4ParameterisationPolyconeZ()
00389 {
00390 }
00391
00392
00393 G4double G4ParameterisationPolyconeZ::GetR(G4double z,
00394 G4double z1, G4double r1,
00395 G4double z2, G4double r2) const
00396 {
00397
00398
00399
00400
00401
00402 return (r1-r2)/(z1-z2)*z + ( r1 - (r1-r2)/(z1-z2)*z1 ) ;
00403 }
00404
00405
00406 G4double G4ParameterisationPolyconeZ::GetRmin(G4double z, G4int nseg) const
00407 {
00408
00409
00410 return GetR(z,
00411 fOrigParamMother->Z_values[nseg],
00412 fOrigParamMother->Rmin[nseg],
00413 fOrigParamMother->Z_values[nseg+1],
00414 fOrigParamMother->Rmin[nseg+1]);
00415 }
00416
00417
00418 G4double G4ParameterisationPolyconeZ::GetRmax(G4double z, G4int nseg) const
00419 {
00420
00421
00422 return GetR(z,
00423 fOrigParamMother->Z_values[nseg],
00424 fOrigParamMother->Rmax[nseg],
00425 fOrigParamMother->Z_values[nseg+1],
00426 fOrigParamMother->Rmax[nseg+1]);
00427 }
00428
00429
00430 G4double G4ParameterisationPolyconeZ::GetMaxParameter() const
00431 {
00432 return std::abs (fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
00433 -fOrigParamMother->Z_values[0]);
00434 }
00435
00436
00437 void G4ParameterisationPolyconeZ::CheckParametersValidity()
00438 {
00439 G4VDivisionParameterisation::CheckParametersValidity();
00440
00441
00442 if( fDivisionType == DivNDIV ) {
00443 if( fnDiv > fOrigParamMother->Num_z_planes-1 ) {
00444 std::ostringstream error;
00445 error << "Configuration not supported." << G4endl
00446 << "Division along Z will be done by splitting in the defined"
00447 << G4endl
00448 << "Z planes, i.e, the number of division would be: "
00449 << fOrigParamMother->Num_z_planes-1
00450 << ", instead of: " << fnDiv << " !";
00451 G4Exception("G4ParameterisationPolyconeZ::CheckParametersValidity()",
00452 "GeomDiv0001", FatalException, error);
00453 }
00454 }
00455
00456
00457
00458 if( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH ) {
00459
00460
00461
00462 G4int isegstart = -1;
00463 G4int isegend = -1;
00464
00465 if ( ! fReflectedSolid ) {
00466
00467 G4double zstart
00468 = fOrigParamMother->Z_values[0] + foffset;
00469 G4double zend
00470 = fOrigParamMother->Z_values[0] + foffset + fnDiv* fwidth;
00471
00472 G4int counter = 0;
00473 while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 ) {
00474
00475 if ( zstart >= fOrigParamMother->Z_values[counter] &&
00476 zstart < fOrigParamMother->Z_values[counter+1] ) {
00477 isegstart = counter;
00478 }
00479
00480 if ( zend > fOrigParamMother->Z_values[counter] &&
00481 zend <= fOrigParamMother->Z_values[counter+1] ) {
00482 isegend = counter;
00483 }
00484 ++counter;
00485 }
00486 }
00487 else {
00488
00489 G4double zstart
00490 = fOrigParamMother->Z_values[0] - foffset;
00491 G4double zend
00492 = fOrigParamMother->Z_values[0] - ( foffset + fnDiv* fwidth);
00493
00494 G4int counter = 0;
00495 while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 ) {
00496
00497 if ( zstart <= fOrigParamMother->Z_values[counter] &&
00498 zstart > fOrigParamMother->Z_values[counter+1] ) {
00499 isegstart = counter;
00500 }
00501
00502 if ( zend < fOrigParamMother->Z_values[counter] &&
00503 zend >= fOrigParamMother->Z_values[counter+1] ) {
00504 isegend = counter;
00505 }
00506 ++counter;
00507 }
00508 }
00509
00510
00511 if ( isegstart != isegend ) {
00512 std::ostringstream message;
00513 message << "Condiguration not supported." << G4endl
00514 << "Division with user defined width." << G4endl
00515 << "Solid " << fmotherSolid->GetName() << G4endl
00516 << "Divided region is not between two z planes.";
00517 G4Exception("G4ParameterisationPolyconeZ::CheckParametersValidity()",
00518 "GeomDiv0001", FatalException, message);
00519 }
00520
00521 fNSegment = isegstart;
00522 }
00523 }
00524
00525
00526 void
00527 G4ParameterisationPolyconeZ::
00528 ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol) const
00529 {
00530 if ( fDivisionType == DivNDIV ) {
00531
00532 G4double posi
00533 = ( fOrigParamMother->Z_values[copyNo]
00534 + fOrigParamMother->Z_values[copyNo+1])/2;
00535 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
00536 }
00537
00538 if ( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH ) {
00539
00540 G4double posi = fOrigParamMother->Z_values[0];
00541
00542 if ( ! fReflectedSolid )
00543 posi += foffset + (2*copyNo + 1) * fwidth/2.;
00544 else
00545 posi -= foffset + (2*copyNo + 1) * fwidth/2.;
00546
00547 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
00548 }
00549
00550
00551
00552 #ifdef G4DIVDEBUG
00553 if( verbose >= 2 )
00554 {
00555 G4cout << " G4ParameterisationPolyconeZ - position: " << posi << G4endl
00556 << " copyNo: " << copyNo << " - foffset: " << foffset/deg
00557 << " - fwidth: " << fwidth/deg << G4endl;
00558 }
00559 #endif
00560
00561 ChangeRotMatrix( physVol );
00562
00563 #ifdef G4DIVDEBUG
00564 if( verbose >= 2 )
00565 {
00566 G4cout << std::setprecision(8) << " G4ParameterisationPolyconeZ "
00567 << copyNo << G4endl
00568 << " Position: " << origin << " - Width: " << fwidth
00569 << " - Axis: " << faxis << G4endl;
00570 }
00571 #endif
00572 }
00573
00574
00575 void
00576 G4ParameterisationPolyconeZ::
00577 ComputeDimensions( G4Polycone& pcone, const G4int copyNo,
00578 const G4VPhysicalVolume* ) const
00579 {
00580
00581
00582 G4PolyconeHistorical origparam;
00583 G4int nz = 2;
00584 origparam.Num_z_planes = nz;
00585 origparam.Start_angle = fOrigParamMother->Start_angle;
00586 origparam.Opening_angle = fOrigParamMother->Opening_angle;
00587
00588
00589 origparam.Z_values = new G4double[nz];
00590 origparam.Rmin = new G4double[nz];
00591 origparam.Rmax = new G4double[nz];
00592
00593 if ( fDivisionType == DivNDIV ) {
00594
00595 G4double posi = (fOrigParamMother->Z_values[copyNo]
00596 + fOrigParamMother->Z_values[copyNo+1])/2;
00597
00598 origparam.Z_values[0] = fOrigParamMother->Z_values[copyNo] - posi;
00599 origparam.Z_values[1] = fOrigParamMother->Z_values[copyNo+1] - posi;
00600 origparam.Rmin[0] = fOrigParamMother->Rmin[copyNo];
00601 origparam.Rmin[1] = fOrigParamMother->Rmin[copyNo+1];
00602 origparam.Rmax[0] = fOrigParamMother->Rmax[copyNo];
00603 origparam.Rmax[1] = fOrigParamMother->Rmax[copyNo+1];
00604 }
00605
00606 if ( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH ) {
00607 if ( ! fReflectedSolid ) {
00608 origparam.Z_values[0] = - fwidth/2.;
00609 origparam.Z_values[1] = fwidth/2.;
00610
00611
00612 G4double posi
00613 = fOrigParamMother->Z_values[0] + foffset + (2*copyNo + 1) * fwidth/2.;
00614
00615
00616 G4double zstart = posi - fwidth/2.;
00617 G4double zend = posi + fwidth/2.;
00618 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
00619 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
00620 origparam.Rmin[1] = GetRmin(zend, fNSegment);
00621 origparam.Rmax[1] = GetRmax(zend, fNSegment);
00622 }
00623 else {
00624 origparam.Z_values[0] = fwidth/2.;
00625 origparam.Z_values[1] = - fwidth/2.;
00626
00627
00628 G4double posi
00629 = fOrigParamMother->Z_values[0] - ( foffset + (2*copyNo + 1) * fwidth/2.);
00630
00631
00632 G4double zstart = posi + fwidth/2.;
00633 G4double zend = posi - fwidth/2.;
00634 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
00635 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
00636 origparam.Rmin[1] = GetRmin(zend, fNSegment);
00637 origparam.Rmax[1] = GetRmax(zend, fNSegment);
00638 }
00639
00640
00641 if ( origparam.Rmin[0] < 0.0 ) origparam.Rmin[0] = 0.0;
00642 if ( origparam.Rmin[nz-1] < 0.0 ) origparam.Rmin[1] = 0.0;
00643 }
00644
00645 pcone.SetOriginalParameters(&origparam);
00646 pcone.Reset();
00647
00648 #ifdef G4DIVDEBUG
00649 if( verbose >= 2 )
00650 {
00651 G4cout << "G4ParameterisationPolyconeZ::ComputeDimensions()" << G4endl
00652 << "-- Parametrised pcone copy-number: " << copyNo << G4endl;
00653 pcone.DumpInfo();
00654 }
00655 #endif
00656 }