#include <G4ParameterisationPolycone.hh>
Inheritance diagram for G4ParameterisationPolyconeZ:
Public Member Functions | |
G4ParameterisationPolyconeZ (EAxis axis, G4int nCopies, G4double offset, G4double step, G4VSolid *motherSolid, DivisionType divType) | |
~G4ParameterisationPolyconeZ () | |
void | CheckParametersValidity () |
G4double | GetMaxParameter () const |
void | ComputeTransformation (const G4int copyNo, G4VPhysicalVolume *physVol) const |
void | ComputeDimensions (G4Polycone &pcone, const G4int copyNo, const G4VPhysicalVolume *physVol) const |
Definition at line 175 of file G4ParameterisationPolycone.hh.
G4ParameterisationPolyconeZ::G4ParameterisationPolyconeZ | ( | EAxis | axis, | |
G4int | nCopies, | |||
G4double | offset, | |||
G4double | step, | |||
G4VSolid * | motherSolid, | |||
DivisionType | divType | |||
) |
Definition at line 352 of file G4ParameterisationPolycone.cc.
References G4VDivisionParameterisation::CalculateNDiv(), CheckParametersValidity(), DivNDIV, DivWIDTH, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4PolyconeHistorical::Num_z_planes, G4VDivisionParameterisation::SetType(), G4VDivisionParameterisation::verbose, and G4PolyconeHistorical::Z_values.
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 }
G4ParameterisationPolyconeZ::~G4ParameterisationPolyconeZ | ( | ) |
void G4ParameterisationPolyconeZ::CheckParametersValidity | ( | ) | [virtual] |
Reimplemented from G4VDivisionParameterisation.
Definition at line 437 of file G4ParameterisationPolycone.cc.
References G4VDivisionParameterisation::CheckParametersValidity(), DivNDIV, DivNDIVandWIDTH, DivWIDTH, FatalException, G4VDivisionParameterisation::fDivisionType, G4VDivisionParameterisation::fmotherSolid, G4VDivisionParameterisation::fnDiv, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fReflectedSolid, G4VDivisionParameterisation::fwidth, G4endl, G4Exception(), G4VSolid::GetName(), G4PolyconeHistorical::Num_z_planes, and G4PolyconeHistorical::Z_values.
Referenced by G4ParameterisationPolyconeZ().
00438 { 00439 G4VDivisionParameterisation::CheckParametersValidity(); 00440 00441 // Division will be following the mother polycone segments 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 // Division will be done within one polycone segment 00457 // with applying given width and offset 00458 if( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH ) { 00459 // Check if divided region does not span over more 00460 // than one z segment 00461 00462 G4int isegstart = -1; // number of the segment containing start position 00463 G4int isegend = -1; // number of the segment containing end position 00464 00465 if ( ! fReflectedSolid ) { 00466 // The start/end position of the divided region 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 // first segment 00475 if ( zstart >= fOrigParamMother->Z_values[counter] && 00476 zstart < fOrigParamMother->Z_values[counter+1] ) { 00477 isegstart = counter; 00478 } 00479 // last segment 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 // The start/end position of the divided region 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 // first segment 00497 if ( zstart <= fOrigParamMother->Z_values[counter] && 00498 zstart > fOrigParamMother->Z_values[counter+1] ) { 00499 isegstart = counter; 00500 } 00501 // last segment 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 }
void G4ParameterisationPolyconeZ::ComputeDimensions | ( | G4Polycone & | pcone, | |
const G4int | copyNo, | |||
const G4VPhysicalVolume * | physVol | |||
) | const [virtual] |
Reimplemented from G4VPVParameterisation.
Definition at line 577 of file G4ParameterisationPolycone.cc.
References DivNDIV, DivNDIVandWIDTH, DivWIDTH, G4VSolid::DumpInfo(), G4VDivisionParameterisation::fDivisionType, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fReflectedSolid, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4PolyconeHistorical::Num_z_planes, G4PolyconeHistorical::Opening_angle, G4Polycone::Reset(), G4PolyconeHistorical::Rmax, G4PolyconeHistorical::Rmin, G4Polycone::SetOriginalParameters(), G4PolyconeHistorical::Start_angle, G4VDivisionParameterisation::verbose, and G4PolyconeHistorical::Z_values.
00579 { 00580 00581 // Define division solid 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 // Define division solid z sections 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 // The position of the centre of copyNo-th mother polycone segment 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 // The position of the centre of copyNo-th division 00612 G4double posi 00613 = fOrigParamMother->Z_values[0] + foffset + (2*copyNo + 1) * fwidth/2.; 00614 00615 // The first and last z sides z values 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 // The position of the centre of copyNo-th division 00628 G4double posi 00629 = fOrigParamMother->Z_values[0] - ( foffset + (2*copyNo + 1) * fwidth/2.); 00630 00631 // The first and last z sides z values 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 // It can happen due to rounding errors 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); // copy values & transfer pointers 00646 pcone.Reset(); // reset to new solid parameters 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 }
void G4ParameterisationPolyconeZ::ComputeTransformation | ( | const G4int | copyNo, | |
G4VPhysicalVolume * | physVol | |||
) | const [virtual] |
Implements G4VDivisionParameterisation.
Definition at line 528 of file G4ParameterisationPolycone.cc.
References G4VDivisionParameterisation::ChangeRotMatrix(), DivNDIV, DivNDIVandWIDTH, DivWIDTH, G4VDivisionParameterisation::faxis, G4VDivisionParameterisation::fDivisionType, G4VDivisionParameterisation::foffset, G4VDivisionParameterisation::fReflectedSolid, G4VDivisionParameterisation::fwidth, G4cout, G4endl, G4VPhysicalVolume::SetTranslation(), G4VDivisionParameterisation::verbose, and G4PolyconeHistorical::Z_values.
00529 { 00530 if ( fDivisionType == DivNDIV ) { 00531 // The position of the centre of copyNo-th mother polycone segment 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 // The position of the centre of copyNo-th division 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 //----- calculate rotation matrix: unit 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 }
G4double G4ParameterisationPolyconeZ::GetMaxParameter | ( | ) | const [virtual] |
Implements G4VDivisionParameterisation.
Definition at line 430 of file G4ParameterisationPolycone.cc.
References G4PolyconeHistorical::Num_z_planes, and G4PolyconeHistorical::Z_values.
00431 { 00432 return std::abs (fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1] 00433 -fOrigParamMother->Z_values[0]); 00434 }