#include <G4tgbVolume.hh>
Public Member Functions | |
G4tgbVolume () | |
~G4tgbVolume () | |
G4tgbVolume (G4tgrVolume *vol) | |
void | ConstructG4Volumes (const G4tgrPlace *place, const G4LogicalVolume *parentLV) |
G4VSolid * | FindOrConstructG4Solid (const G4tgrSolid *vol) |
G4LogicalVolume * | ConstructG4LogVol (const G4VSolid *solid) |
G4VPhysicalVolume * | ConstructG4PhysVol (const G4tgrPlace *place, const G4LogicalVolume *currentLV, const G4LogicalVolume *parentLV) |
void | SetCutsInRange (G4LogicalVolume *logvol, std::map< G4String, G4double > cuts) |
void | SetCutsInEnergy (G4LogicalVolume *logvol, std::map< G4String, G4double > cuts) |
void | CheckNoSolidParams (const G4String &solidType, const unsigned int NoParamExpected, const unsigned int NoParam) |
G4VSolid * | BuildSolidForDivision (G4VSolid *parentSolid, EAxis axis) |
const G4String & | GetName () const |
G4bool | GetVisibility () const |
const G4double * | GetColour () const |
Definition at line 66 of file G4tgbVolume.hh.
G4tgbVolume::G4tgbVolume | ( | ) |
G4tgbVolume::~G4tgbVolume | ( | ) |
G4tgbVolume::G4tgbVolume | ( | G4tgrVolume * | vol | ) |
Definition at line 1306 of file G4tgbVolume.cc.
References FatalException, G4cout, G4endl, G4Exception(), G4VSolid::GetEntityType(), G4VSolid::GetExtent(), G4GeometryTolerance::GetInstance(), GetName(), G4tgrMessenger::GetVerboseLevel(), G4VisExtent::GetXmax(), G4VisExtent::GetXmin(), G4VisExtent::GetYmax(), G4VisExtent::GetYmin(), G4VisExtent::GetZmax(), G4VisExtent::GetZmin(), kXAxis, G4PolyhedraHistorical::Num_z_planes, G4PolyconeHistorical::Num_z_planes, G4PolyhedraHistorical::Rmax, G4PolyconeHistorical::Rmax, G4PolyhedraHistorical::Rmin, G4PolyconeHistorical::Rmin, G4PolyhedraHistorical::Z_values, and G4PolyconeHistorical::Z_values.
Referenced by ConstructG4PhysVol().
01307 { 01308 G4VSolid* solid=0; 01309 G4double redf = (parentSolid->GetExtent().GetXmax()-parentSolid->GetExtent().GetXmin()); 01310 redf = std::min(redf,parentSolid->GetExtent().GetYmax()-parentSolid->GetExtent().GetYmin()); 01311 redf = std::min(redf,parentSolid->GetExtent().GetZmax()-parentSolid->GetExtent().GetZmin()); 01312 redf *= 0.001; //make daugther much smaller, to fit in parent 01313 01314 if( parentSolid->GetEntityType() == "G4Box" ) 01315 { 01316 G4Box* psolid = (G4Box*)(parentSolid); 01317 solid = new G4Box(GetName(), psolid->GetXHalfLength()*redf, 01318 psolid->GetZHalfLength()*redf, 01319 psolid->GetZHalfLength()*redf); 01320 } 01321 else if ( parentSolid->GetEntityType() == "G4Tubs" ) 01322 { 01323 G4Tubs* psolid = (G4Tubs*)(parentSolid); 01324 solid = new G4Tubs( GetName(), psolid->GetInnerRadius()*redf, 01325 psolid->GetOuterRadius()*redf, 01326 psolid->GetZHalfLength()*redf, 01327 psolid->GetSPhi(), psolid->GetDPhi()); 01328 } 01329 else if ( parentSolid->GetEntityType() == "G4Cons" ) 01330 { 01331 G4Cons* psolid = (G4Cons*)(parentSolid); 01332 solid = new G4Cons( GetName(), psolid->GetInnerRadiusMinusZ()*redf, 01333 psolid->GetOuterRadiusMinusZ()*redf, 01334 psolid->GetInnerRadiusPlusZ()*redf, 01335 psolid->GetOuterRadiusPlusZ()*redf, 01336 psolid->GetZHalfLength()*redf, 01337 psolid->GetSPhi(), psolid->GetDPhi()); 01338 } 01339 else if ( parentSolid->GetEntityType() == "G4Trd" ) 01340 { 01341 G4Trd* psolid = (G4Trd*)(parentSolid); 01342 G4double mpDx1 = psolid->GetXHalfLength1(); 01343 G4double mpDx2 = psolid->GetXHalfLength2(); 01344 01345 if( axis == kXAxis && std::fabs(mpDx1 - mpDx2) > G4GeometryTolerance::GetInstance()->GetSurfaceTolerance() ) 01346 { 01347 solid = new G4Trap( GetName(), psolid->GetZHalfLength()*redf, 01348 psolid->GetYHalfLength1()*redf, 01349 psolid->GetXHalfLength2()*redf, 01350 psolid->GetXHalfLength1()*redf ); 01351 } 01352 else 01353 { 01354 solid = new G4Trd( GetName(), psolid->GetXHalfLength1()*redf, 01355 psolid->GetXHalfLength2()*redf, 01356 psolid->GetYHalfLength1()*redf, 01357 psolid->GetYHalfLength2()*redf, 01358 psolid->GetZHalfLength()*redf); 01359 } 01360 01361 } 01362 else if ( parentSolid->GetEntityType() == "G4Para" ) 01363 { 01364 G4Para* psolid = (G4Para*)(parentSolid); 01365 solid = new G4Para( GetName(), psolid->GetXHalfLength()*redf, 01366 psolid->GetYHalfLength()*redf, 01367 psolid->GetZHalfLength()*redf, 01368 std::atan(psolid->GetTanAlpha()), 01369 psolid->GetSymAxis().theta(), 01370 psolid->GetSymAxis().phi() ); 01371 } 01372 else if ( parentSolid->GetEntityType() == "G4Polycone" ) 01373 { 01374 G4Polycone* psolid = (G4Polycone*)(parentSolid); 01375 G4PolyconeHistorical origParam = *(psolid->GetOriginalParameters()); 01376 for( G4int ii = 0; ii < origParam.Num_z_planes; ii++ ) 01377 { 01378 origParam.Rmin[ii] = origParam.Rmin[ii]*redf; 01379 origParam.Rmax[ii] = origParam.Rmax[ii]*redf; 01380 } 01381 solid = new G4Polycone( GetName(), psolid->GetStartPhi(), 01382 psolid->GetEndPhi(), 01383 origParam.Num_z_planes, origParam.Z_values, 01384 origParam.Rmin, origParam.Rmax); 01385 01386 } 01387 else if ( parentSolid->GetEntityType() == "G4Polyhedra" ) 01388 { 01389 G4Polyhedra* psolid = (G4Polyhedra*)(parentSolid); 01390 G4PolyhedraHistorical origParam = *(psolid->GetOriginalParameters()); 01391 for( G4int ii = 0; ii < origParam.Num_z_planes; ii++ ) 01392 { 01393 origParam.Rmin[ii] = origParam.Rmin[ii]*redf; 01394 origParam.Rmax[ii] = origParam.Rmax[ii]*redf; 01395 } 01396 solid = new G4Polyhedra( GetName(), psolid->GetStartPhi(), 01397 psolid->GetEndPhi(), 01398 psolid->GetNumSide(), 01399 origParam.Num_z_planes, origParam.Z_values, 01400 origParam.Rmin, origParam.Rmax); 01401 } 01402 else 01403 { 01404 G4String ErrMessage = "Solid type not supported. VOLUME= " + GetName() 01405 + " Solid type= " + parentSolid->GetEntityType() + "\n" 01406 + "Only supported types are: G4Box, G4Tubs, G4Cons," 01407 + " G4Trd, G4Para, G4Polycone, G4Polyhedra."; 01408 G4Exception("G4tgbVolume::BuildSolidForDivision()", "NotImplemented", 01409 FatalException, ErrMessage); 01410 return 0; 01411 } 01412 01413 #ifdef G4VERBOSE 01414 if( G4tgrMessenger::GetVerboseLevel() >= 1 ) 01415 { 01416 G4cout << " Constructing new G4Solid for division: " 01417 << *solid << G4endl; 01418 } 01419 #endif 01420 return solid; 01421 }
void G4tgbVolume::CheckNoSolidParams | ( | const G4String & | solidType, | |
const unsigned int | NoParamExpected, | |||
const unsigned int | NoParam | |||
) |
Definition at line 828 of file G4tgbVolume.cc.
References G4UIcommand::ConvertToString(), FatalException, and G4Exception().
Referenced by FindOrConstructG4Solid().
00831 { 00832 if( NoParamExpected != NoParam ) 00833 { 00834 G4String Err1 = "Solid type " + solidType + " should have "; 00835 G4String Err2 = G4UIcommand::ConvertToString(G4int(NoParamExpected)) 00836 + " parameters,\n"; 00837 G4String Err3 = "and it has " 00838 + G4UIcommand::ConvertToString(G4int(NoParam)); 00839 G4String ErrMessage = Err1 + Err2 + Err3 + " !"; 00840 G4Exception("G4tgbVolume::CheckNoSolidParams()", "InvalidSetup", 00841 FatalException, ErrMessage); 00842 } 00843 }
G4LogicalVolume * G4tgbVolume::ConstructG4LogVol | ( | const G4VSolid * | solid | ) |
Definition at line 847 of file G4tgbVolume.cc.
References FatalException, G4tgbMaterialMgr::FindOrBuildG4Material(), G4cout, G4endl, G4Exception(), GetColour(), G4tgbMaterialMgr::GetInstance(), G4tgrVolume::GetMaterialName(), G4LogicalVolume::GetName(), G4Material::GetName(), GetName(), G4tgrMessenger::GetVerboseLevel(), GetVisibility(), G4VisAttributes::SetColour(), G4LogicalVolume::SetVisAttributes(), and G4VisAttributes::SetVisibility().
Referenced by ConstructG4Volumes().
00848 { 00849 G4LogicalVolume* logvol; 00850 00851 #ifdef G4VERBOSE 00852 if( G4tgrMessenger::GetVerboseLevel() >= 2 ) 00853 { 00854 G4cout << " G4tgbVolume::ConstructG4LogVol() - " << GetName() << G4endl; 00855 } 00856 #endif 00857 00858 //----------- Get the material first 00859 G4Material* mate = G4tgbMaterialMgr::GetInstance() 00860 ->FindOrBuildG4Material( theTgrVolume->GetMaterialName() ); 00861 if( mate == 0 ) 00862 { 00863 G4String ErrMessage = "Material not found " 00864 + theTgrVolume->GetMaterialName() 00865 + " for volume " + GetName() + "."; 00866 G4Exception("G4tgbVolume::ConstructG4LogVol()", "InvalidSetup", 00867 FatalException, ErrMessage); 00868 } 00869 #ifdef G4VERBOSE 00870 if( G4tgrMessenger::GetVerboseLevel() >= 2 ) 00871 { 00872 G4cout << " G4tgbVolume::ConstructG4LogVol() -" 00873 << " Material constructed: " << mate->GetName() << G4endl; 00874 } 00875 #endif 00876 00877 //---------- Construct the LV 00878 logvol = new G4LogicalVolume( const_cast<G4VSolid*>(solid), 00879 const_cast<G4Material*>(mate), GetName() ); 00880 00881 #ifdef G4VERBOSE 00882 if( G4tgrMessenger::GetVerboseLevel() >= 1 ) 00883 { 00884 G4cout << " Constructing new G4LogicalVolume: " 00885 << logvol->GetName() << " mate " << mate->GetName() << G4endl; 00886 } 00887 #endif 00888 00889 //---------- Set visibility and colour 00890 if( !GetVisibility() || GetColour()[0] != -1 ) 00891 { 00892 G4VisAttributes* visAtt = new G4VisAttributes(); 00893 #ifdef G4VERBOSE 00894 if( G4tgrMessenger::GetVerboseLevel() >= 1 ) 00895 { 00896 G4cout << " Constructing new G4VisAttributes: " 00897 << *visAtt << G4endl; 00898 } 00899 #endif 00900 00901 if( !GetVisibility() ) 00902 { 00903 visAtt->SetVisibility( GetVisibility() ); 00904 } 00905 else if( GetColour()[0] != -1 ) 00906 { 00907 // this else should not be necessary, because if the visibility 00908 // is set to off, colour should have no effect. But it does not 00909 // work: if you set colour and vis off, it is visualized!?!?!? 00910 00911 const G4double* col = GetColour(); 00912 if( col[3] == -1. ) 00913 { 00914 visAtt->SetColour( G4Colour(col[0],col[1],col[2])); 00915 } 00916 else 00917 { 00918 visAtt->SetColour( G4Colour(col[0],col[1],col[2],col[3])); 00919 } 00920 } 00921 logvol->SetVisAttributes(visAtt); 00922 } 00923 00924 #ifdef G4VERBOSE 00925 if( G4tgrMessenger::GetVerboseLevel() >= 2 ) 00926 { 00927 G4cout << " G4tgbVolume::ConstructG4LogVol() -" 00928 << " Created logical volume: " << GetName() << G4endl; 00929 } 00930 #endif 00931 00932 return logvol; 00933 }
G4VPhysicalVolume * G4tgbVolume::ConstructG4PhysVol | ( | const G4tgrPlace * | place, | |
const G4LogicalVolume * | currentLV, | |||
const G4LogicalVolume * | parentLV | |||
) |
Definition at line 938 of file G4tgbVolume.cc.
References G4AssemblyVolume::AddPlacedVolume(), BuildSolidForDivision(), ConstructG4Volumes(), DivByNdiv, DivByNdivAndWidth, DivByWidth, FatalException, G4tgbVolumeMgr::FindG4LogVol(), G4tgbMaterialMgr::FindOrBuildG4Material(), G4tgbRotationMatrixMgr::FindOrBuildG4RotMatrix(), G4tgbVolumeMgr::FindVolume(), G4cerr, G4cout, G4endl, G4Exception(), G4tgrPlaceDivRep::GetAxis(), G4tgbPlaceParameterisation::GetAxis(), G4tgrVolume::GetCheckOverlaps(), G4tgrPlace::GetCopyNo(), G4tgrPlaceDivRep::GetDivType(), G4tgbVolumeMgr::GetInstance(), G4tgbMaterialMgr::GetInstance(), G4tgbRotationMatrixMgr::GetInstance(), G4tgrVolume::GetMaterialName(), G4LogicalVolume::GetName(), G4VPhysicalVolume::GetName(), GetName(), G4tgbPlaceParameterisation::GetNCopies(), G4tgrPlaceDivRep::GetNDiv(), G4tgrPlaceDivRep::GetOffset(), G4tgrPlace::GetPlacement(), G4tgrPlaceSimple::GetRotMatName(), G4LogicalVolume::GetSolid(), G4tgrPlace::GetType(), G4tgrVolume::GetType(), G4tgrMessenger::GetVerboseLevel(), G4tgrPlaceDivRep::GetWidth(), G4ReflectionFactory::Instance(), G4AssemblyVolume::MakeImprint(), and G4ReflectionFactory::Place().
Referenced by ConstructG4Volumes().
00941 { 00942 G4VPhysicalVolume* physvol = 0; 00943 G4int copyNo; 00944 00945 //----- Case of placement of top volume 00946 if( place == 0 ) 00947 { 00948 #ifdef G4VERBOSE 00949 if( G4tgrMessenger::GetVerboseLevel() >= 2 ) 00950 { 00951 G4cout << " G4tgbVolume::ConstructG4PhysVol() - World: " 00952 << GetName() << G4endl; 00953 } 00954 #endif 00955 physvol = new G4PVPlacement(0, G4ThreeVector(), 00956 const_cast<G4LogicalVolume*>(currentLV), 00957 GetName(), 0, false, 0, 00958 theTgrVolume->GetCheckOverlaps()); 00959 #ifdef G4VERBOSE 00960 if( G4tgrMessenger::GetVerboseLevel() >= 1 ) 00961 { 00962 G4cout << " Constructing new : G4PVPlacement " 00963 << physvol->GetName() << G4endl; 00964 } 00965 #endif 00966 } 00967 else 00968 { 00969 copyNo = place->GetCopyNo(); 00970 00971 #ifdef G4VERBOSE 00972 if( G4tgrMessenger::GetVerboseLevel() >= 2 ) 00973 { 00974 G4cout << " G4tgbVolume::ConstructG4PhysVol() - " << GetName() << G4endl 00975 << " inside " << parentLV->GetName() << " copy No: " << copyNo 00976 << " of type= " << theTgrVolume->GetType() << G4endl 00977 << " placement type= " << place->GetType() << G4endl; 00978 } 00979 #endif 00980 00981 if( theTgrVolume->GetType() == "VOLSimple" ) 00982 { 00983 //----- Get placement 00984 #ifdef G4VERBOSE 00985 if( G4tgrMessenger::GetVerboseLevel() >= 2 ) 00986 { 00987 G4cout << " G4tgbVolume::ConstructG4PhysVol() - Placement type = " 00988 << place->GetType() << G4endl; 00989 } 00990 #endif 00991 00992 //--------------- If it is G4tgrPlaceSimple 00993 if( place->GetType() == "PlaceSimple" ) 00994 { 00995 //----- Get rotation matrix 00996 G4tgrPlaceSimple* placeSimple = (G4tgrPlaceSimple*)place; 00997 G4String rmName = placeSimple->GetRotMatName(); 00998 00999 G4RotationMatrix* rotmat = G4tgbRotationMatrixMgr::GetInstance() 01000 ->FindOrBuildG4RotMatrix( rmName ); 01001 //----- Place volume in mother 01002 G4double check = (rotmat->colX().cross(rotmat->colY()))*rotmat->colZ(); 01003 G4double tol = 1.0e-3; 01004 //---- Check that matrix is ortogonal 01005 if (1-std::abs(check)>tol) 01006 { 01007 G4cerr << " Matrix : " << rmName << " " << rotmat->colX() 01008 << " " << rotmat->colY() << " " << rotmat->colZ() << G4endl 01009 << " product x X y * z = " << check << " x X y " 01010 << rotmat->colX().cross(rotmat->colY()) << G4endl; 01011 G4String ErrMessage = "Rotation is not ortogonal " + rmName + " !"; 01012 G4Exception("G4tgbVolume::ConstructG4PhysVol()", 01013 "InvalidSetup", FatalException, ErrMessage); 01014 //---- Check if it is reflection 01015 } 01016 else if (1+check<=tol) 01017 { 01018 G4Translate3D transl = place->GetPlacement(); 01019 G4Transform3D trfrm = transl * G4Rotate3D(*rotmat); 01020 physvol = (G4ReflectionFactory::Instance()->Place(trfrm, GetName(), 01021 const_cast<G4LogicalVolume*>(currentLV), 01022 const_cast<G4LogicalVolume*>(parentLV), 01023 false, copyNo, false )).first; 01024 } 01025 else 01026 { 01027 #ifdef G4VERBOSE 01028 if( G4tgrMessenger::GetVerboseLevel() >= 1 ) 01029 { 01030 G4cout << "Construction new G4VPhysicalVolume" 01031 << " through G4ReflectionFactory " << GetName() 01032 << " in volume " << parentLV->GetName() 01033 << " copyNo " << copyNo 01034 << " position " << place->GetPlacement() 01035 << " ROT " << rotmat->colX() 01036 << " " << rotmat->colY() 01037 << " " << rotmat->colZ() << G4endl; 01038 } 01039 #endif 01040 physvol = new G4PVPlacement( rotmat, place->GetPlacement(), 01041 const_cast<G4LogicalVolume*>(currentLV), 01042 GetName(), 01043 const_cast<G4LogicalVolume*>(parentLV), 01044 false, copyNo, 01045 theTgrVolume->GetCheckOverlaps()); 01046 } 01047 01048 //--------------- If it is G4tgrPlaceParam 01049 } 01050 else if( place->GetType() == "PlaceParam" ) 01051 { 01052 G4tgrPlaceParameterisation* dp = (G4tgrPlaceParameterisation*)(place); 01053 01054 //----- See what parameterisation type 01055 #ifdef G4VERBOSE 01056 if( G4tgrMessenger::GetVerboseLevel() >= 2 ) 01057 { 01058 G4cout << " G4tgbVolume::ConstructG4PhysVol() -" << G4endl 01059 << " param: " << GetName() << " in " << parentLV->GetName() 01060 << " param type= " << dp->GetParamType() << G4endl; 01061 } 01062 #endif 01063 01064 G4tgbPlaceParameterisation * param=0; 01065 01066 if( (dp->GetParamType() == "CIRCLE") 01067 || (dp->GetParamType() == "CIRCLE_XY") 01068 || (dp->GetParamType() == "CIRCLE_XZ") 01069 || (dp->GetParamType() == "CIRCLE_YZ") ) 01070 { 01071 param = new G4tgbPlaceParamCircle(dp); 01072 01073 } 01074 else if( (dp->GetParamType() == "LINEAR") 01075 || (dp->GetParamType() == "LINEAR_X") 01076 || (dp->GetParamType() == "LINEAR_Y") 01077 || (dp->GetParamType() == "LINEAR_Z") ) 01078 { 01079 param = new G4tgbPlaceParamLinear(dp); 01080 01081 } 01082 else if( (dp->GetParamType() == "SQUARE") 01083 || (dp->GetParamType() == "SQUARE_XY") 01084 || (dp->GetParamType() == "SQUARE_XZ") 01085 || (dp->GetParamType() == "SQUARE_YZ") ) 01086 { 01087 param = new G4tgbPlaceParamSquare(dp); 01088 } 01089 else 01090 { 01091 G4String ErrMessage = "Parameterisation has wrong type, TYPE: " 01092 + G4String(dp->GetParamType()) + " !"; 01093 G4Exception("G4tgbVolume::ConstructG4PhysVol", "WrongArgument", 01094 FatalException, ErrMessage); 01095 return 0; 01096 } 01097 #ifdef G4VERBOSE 01098 if( G4tgrMessenger::GetVerboseLevel() >= 1 ) 01099 { 01100 G4cout << " G4tgbVolume::ConstructG4PhysVol() -" << G4endl 01101 << " New G4PVParameterised: " << GetName() << " vol " 01102 << currentLV->GetName() << " in vol " << parentLV->GetName() 01103 << " axis " << param->GetAxis() << " nCopies " 01104 << param->GetNCopies() << G4endl; 01105 } 01106 #endif 01107 physvol = new G4PVParameterised(GetName(), 01108 const_cast<G4LogicalVolume*>(currentLV), 01109 const_cast<G4LogicalVolume*>(parentLV), 01110 EAxis(param->GetAxis()), 01111 param->GetNCopies(), param); 01112 #ifdef G4VERBOSE 01113 if( G4tgrMessenger::GetVerboseLevel() >= 1 ) 01114 { 01115 G4cout << " Constructing new G4PVParameterised: " 01116 << physvol->GetName() << " in volume " << parentLV->GetName() 01117 << " N copies " << param->GetNCopies() 01118 << " axis " << param->GetAxis() << G4endl; 01119 } 01120 #endif 01121 01122 } 01123 else if( place->GetType() == "PlaceReplica" ) 01124 { 01125 //--------------- If it is PlaceReplica 01126 G4tgrPlaceDivRep* dpr = (G4tgrPlaceDivRep*)place; 01127 01128 #ifdef G4VERBOSE 01129 if( G4tgrMessenger::GetVerboseLevel() >= 2 ) 01130 { 01131 G4cout << " G4tgbVolume::ConstructG4PhysVol() -" << G4endl 01132 << " replica" << " " << currentLV->GetName() 01133 << " in " << parentLV->GetName() 01134 << " NDiv " << dpr->GetNDiv() << " Width " << dpr->GetWidth() 01135 << " offset " << dpr->GetOffset() << G4endl; 01136 } 01137 #endif 01138 physvol = new G4PVReplica(GetName(), 01139 const_cast<G4LogicalVolume*>(currentLV), 01140 const_cast<G4LogicalVolume*>(parentLV), 01141 EAxis(dpr->GetAxis()), dpr->GetNDiv(), 01142 dpr->GetWidth(), dpr->GetOffset()); 01143 #ifdef G4VERBOSE 01144 if( G4tgrMessenger::GetVerboseLevel() >= 1 ) 01145 { 01146 G4cout << " Constructing new G4PVReplica: " 01147 << currentLV->GetName() 01148 << " in " << parentLV->GetName() 01149 << " NDiv " << dpr->GetNDiv() << " Width " << dpr->GetWidth() 01150 << " offset " << dpr->GetOffset() << G4endl; 01151 } 01152 #endif 01153 } 01154 } 01155 else if( theTgrVolume->GetType() == "VOLDivision" ) 01156 { 01157 G4tgrVolumeDivision* volr = (G4tgrVolumeDivision*)theTgrVolume; 01158 G4tgrPlaceDivRep* placeDiv = volr->GetPlaceDivision() ; 01159 G4VSolid* solid = BuildSolidForDivision( parentLV->GetSolid(), placeDiv->GetAxis() ); 01160 G4Material* mate = G4tgbMaterialMgr::GetInstance() 01161 ->FindOrBuildG4Material( theTgrVolume->GetMaterialName() ); 01162 G4LogicalVolume* divLV = new G4LogicalVolume(solid, 01163 const_cast<G4Material*>(mate), 01164 GetName() ); 01165 #ifdef G4VERBOSE 01166 if( G4tgrMessenger::GetVerboseLevel() >= 1 ) 01167 { 01168 G4cout << " Constructed new G4LogicalVolume for division: " 01169 << divLV->GetName() << " mate " << mate->GetName() << G4endl; 01170 } 01171 #endif 01172 01173 G4DivType divType = placeDiv->GetDivType(); 01174 switch (divType) 01175 { 01176 case DivByNdiv: 01177 physvol = new G4PVDivision(GetName(), (G4LogicalVolume*)divLV, 01178 const_cast<G4LogicalVolume*>(parentLV), 01179 placeDiv->GetAxis(), placeDiv->GetNDiv(), 01180 placeDiv->GetOffset()); 01181 #ifdef G4VERBOSE 01182 if( G4tgrMessenger::GetVerboseLevel() >= 1 ) 01183 { 01184 G4cout << " Constructing new G4PVDivision by number of divisions: " 01185 << GetName() << " in " << parentLV->GetName() 01186 << " axis " << placeDiv->GetAxis() 01187 << " Ndiv " << placeDiv->GetNDiv() 01188 << " offset " << placeDiv->GetOffset() << G4endl; 01189 } 01190 #endif 01191 break; 01192 case DivByWidth: 01193 physvol = new G4PVDivision(GetName(), (G4LogicalVolume*)divLV, 01194 const_cast<G4LogicalVolume*>(parentLV), 01195 placeDiv->GetAxis(), placeDiv->GetWidth(), 01196 placeDiv->GetOffset()); 01197 #ifdef G4VERBOSE 01198 if( G4tgrMessenger::GetVerboseLevel() >= 1 ) 01199 { 01200 G4cout << " Constructing new G4PVDivision by width: " 01201 << GetName() << " in " << parentLV->GetName() 01202 << " axis " << placeDiv->GetAxis() 01203 << " width " << placeDiv->GetWidth() 01204 << " offset " << placeDiv->GetOffset() << G4endl; 01205 } 01206 #endif 01207 break; 01208 case DivByNdivAndWidth: 01209 physvol = new G4PVDivision(GetName(), (G4LogicalVolume*)divLV, 01210 const_cast<G4LogicalVolume*>(parentLV), 01211 placeDiv->GetAxis(), placeDiv->GetNDiv(), 01212 placeDiv->GetWidth(), 01213 placeDiv->GetOffset()); 01214 #ifdef G4VERBOSE 01215 if( G4tgrMessenger::GetVerboseLevel() >= 1 ) 01216 { 01217 G4cout << " Constructing new G4PVDivision by width" 01218 << " and number of divisions: " 01219 << GetName() << " in " << parentLV->GetName() 01220 << " axis " << placeDiv->GetAxis() 01221 << " Ndiv " << placeDiv->GetNDiv() 01222 << " width " << placeDiv->GetWidth() 01223 << " offset " << placeDiv->GetOffset() << G4endl; 01224 } 01225 #endif 01226 break; 01227 } 01228 } 01229 else if( theTgrVolume->GetType() == "VOLAssembly" ) 01230 { 01231 // Define one layer as one assembly volume 01232 G4tgrVolumeAssembly * tgrAssembly = (G4tgrVolumeAssembly *)theTgrVolume; 01233 01234 if( !theG4AssemblyVolume ) 01235 { 01236 theG4AssemblyVolume = new G4AssemblyVolume; 01237 #ifdef G4VERBOSE 01238 if( G4tgrMessenger::GetVerboseLevel() >= 1 ) 01239 { 01240 G4cout << " Constructing new G4AssemblyVolume: " 01241 << " number of assembly components " 01242 << tgrAssembly->GetNoComponents() << G4endl; 01243 } 01244 #endif 01245 G4tgbVolumeMgr* g4vmgr = G4tgbVolumeMgr::GetInstance(); 01246 for( G4int ii = 0; ii < tgrAssembly->GetNoComponents(); ii++ ) 01247 { 01248 // Rotation and translation of a plate inside the assembly 01249 01250 G4ThreeVector transl = tgrAssembly->GetComponentPos(ii); 01251 G4String rmName = tgrAssembly->GetComponentRM(ii); 01252 G4RotationMatrix* rotmat = G4tgbRotationMatrixMgr::GetInstance() 01253 ->FindOrBuildG4RotMatrix( rmName ); 01254 01255 //----- Get G4LogicalVolume of component 01256 G4String lvname = tgrAssembly->GetComponentName(ii); 01257 G4LogicalVolume* logvol = g4vmgr->FindG4LogVol( lvname); 01258 if( logvol == 0 ) 01259 { 01260 g4vmgr->FindVolume( lvname )->ConstructG4Volumes( 0, 0); 01261 logvol = g4vmgr->FindG4LogVol( lvname, true ); 01262 } 01263 // Fill the assembly by the plates 01264 theG4AssemblyVolume->AddPlacedVolume( logvol, transl, rotmat ); 01265 #ifdef G4VERBOSE 01266 if( G4tgrMessenger::GetVerboseLevel() >= 1 ) 01267 { 01268 G4cout << " G4AssemblyVolume->AddPlacedVolume " << ii 01269 << " " << logvol->GetName() 01270 << " translation " << transl 01271 << " rotmat " << rotmat->colX() 01272 << " " << rotmat->colY() 01273 << " " << rotmat->colZ() << G4endl; 01274 } 01275 #endif 01276 } 01277 } 01278 01279 // Rotation and Translation of the assembly inside the world 01280 01281 G4tgrPlaceSimple* placeSimple = (G4tgrPlaceSimple*)place; 01282 G4String rmName = placeSimple->GetRotMatName(); 01283 G4RotationMatrix* rotmat = G4tgbRotationMatrixMgr::GetInstance() 01284 ->FindOrBuildG4RotMatrix( rmName ); 01285 G4ThreeVector transl = place->GetPlacement(); 01286 01287 G4LogicalVolume* parentLV_nonconst = 01288 const_cast<G4LogicalVolume*>(parentLV); 01289 theG4AssemblyVolume->MakeImprint( parentLV_nonconst, transl, rotmat ); 01290 01291 } 01292 else // If it is G4tgrVolumeAssembly 01293 { 01294 G4String ErrMessage = "Volume type not supported: " 01295 + theTgrVolume->GetType() + ", sorry..."; 01296 G4Exception("G4tgbVolume::ConstructG4PhysVol()", "NotImplemented", 01297 FatalException, ErrMessage); 01298 } 01299 } 01300 01301 return physvol; 01302 }
void G4tgbVolume::ConstructG4Volumes | ( | const G4tgrPlace * | place, | |
const G4LogicalVolume * | parentLV | |||
) |
Definition at line 137 of file G4tgbVolume.cc.
References ConstructG4LogVol(), ConstructG4PhysVol(), G4tgbVolumeMgr::FindG4LogVol(), FindOrConstructG4Solid(), G4tgbVolumeMgr::FindVolume(), G4cout, G4endl, G4tgrVolumeMgr::GetChildren(), G4tgrPlace::GetCopyNo(), G4tgrVolumeMgr::GetInstance(), G4tgbVolumeMgr::GetInstance(), G4VPhysicalVolume::GetLogicalVolume(), G4tgrVolume::GetName(), G4LogicalVolume::GetName(), GetName(), G4tgrVolume::GetSolid(), G4tgrVolume::GetType(), G4tgrMessenger::GetVerboseLevel(), G4tgrPlace::GetVolume(), G4tgbVolumeMgr::RegisterChildParentLVs(), and G4tgbVolumeMgr::RegisterMe().
Referenced by G4tgbDetectorConstruction::Construct(), G4tgbDetectorBuilder::ConstructDetector(), and ConstructG4PhysVol().
00139 { 00140 #ifdef G4VERBOSE 00141 if( G4tgrMessenger::GetVerboseLevel() >= 2 ) 00142 { 00143 G4cout << G4endl << "@@@ G4tgbVolume::ConstructG4Volumes - " << GetName() << G4endl; 00144 if( place && parentLV ) G4cout << " place in LV " << parentLV->GetName() << G4endl; 00145 } 00146 #endif 00147 G4tgbVolumeMgr* g4vmgr = G4tgbVolumeMgr::GetInstance(); 00148 G4LogicalVolume* logvol = g4vmgr->FindG4LogVol( GetName() ); 00149 G4bool bFirstCopy = false; 00150 if( logvol == 0 ) 00151 { 00152 bFirstCopy = true; 00153 if( theTgrVolume->GetType() != "VOLDivision" ) 00154 { 00155 //--- If first time build solid and LogVol 00156 G4VSolid* solid = FindOrConstructG4Solid( theTgrVolume->GetSolid() ); 00157 if( solid != 0 ) // for G4AssemblyVolume it is 0 00158 { 00159 g4vmgr->RegisterMe( solid ); 00160 logvol = ConstructG4LogVol( solid ); 00161 g4vmgr->RegisterMe( logvol ); 00162 g4vmgr->RegisterChildParentLVs( logvol, parentLV ); 00163 } 00164 } 00165 else 00166 { 00167 return; 00168 } 00169 } 00170 //--- Construct PhysVol 00171 G4VPhysicalVolume* physvol = ConstructG4PhysVol( place, logvol, parentLV ); 00172 if( physvol != 0 ) // 0 for G4AssemblyVolumes 00173 { 00174 g4vmgr->RegisterMe( physvol ); 00175 00176 if( logvol == 0 ) // case of divisions 00177 { 00178 logvol = physvol->GetLogicalVolume(); 00179 } 00180 } 00181 else 00182 { 00183 return; 00184 } 00185 00186 //--- If first copy build children placements in this LogVol 00187 if(bFirstCopy) 00188 { 00189 std::pair<G4mmapspl::iterator, G4mmapspl::iterator> children 00190 = G4tgrVolumeMgr::GetInstance()->GetChildren( GetName() ); 00191 G4mmapspl::iterator cite; 00192 for( cite = children.first; cite != children.second; cite++ ) 00193 { 00194 //----- Call G4tgrPlace ->constructG4Volumes 00195 //---- find G4tgbVolume corresponding to the G4tgrVolume 00196 // pointed by G4tgrPlace 00197 G4tgrPlace* pl = const_cast<G4tgrPlace*>((*cite).second); 00198 G4tgbVolume* svol = g4vmgr->FindVolume( pl->GetVolume()->GetName() ); 00199 //--- find copyNo 00200 #ifdef G4VERBOSE 00201 if( G4tgrMessenger::GetVerboseLevel() >= 2 ) 00202 { 00203 G4cout << " G4tgbVolume::ConstructG4Volumes - construct daughter " << pl->GetVolume()->GetName() << " # " << pl->GetCopyNo() << G4endl; 00204 } 00205 #endif 00206 svol->ConstructG4Volumes( pl, logvol ); 00207 } 00208 } 00209 00210 }
G4VSolid * G4tgbVolume::FindOrConstructG4Solid | ( | const G4tgrSolid * | vol | ) |
Definition at line 215 of file G4tgbVolume.cc.
References ABSOLUTE, G4TessellatedSolid::AddFacet(), CheckNoSolidParams(), G4UIcommand::ConvertToString(), FatalException, G4tgbVolumeMgr::FindG4Solid(), G4tgbRotationMatrixMgr::FindOrBuildG4RotMatrix(), G4cout, G4endl, G4Exception(), G4GeometryTolerance::GetAngularTolerance(), G4VSolid::GetEntityType(), G4tgbRotationMatrixMgr::GetInstance(), G4tgbVolumeMgr::GetInstance(), G4GeometryTolerance::GetInstance(), G4tgrSolid::GetName(), G4tgrSolid::GetSolidParams(), G4tgrSolid::GetType(), G4tgrMessenger::GetVerboseLevel(), G4INCL::Math::pi, and RELATIVE.
Referenced by ConstructG4Volumes().
00216 { 00217 G4double angularTolerance = G4GeometryTolerance::GetInstance() 00218 ->GetAngularTolerance(); 00219 00220 if( sol == 0 ) { return 0; } 00221 00222 #ifdef G4VERBOSE 00223 if( G4tgrMessenger::GetVerboseLevel() >= 2 ) 00224 { 00225 G4cout << " G4tgbVolume::FindOrConstructG4Solid():" << G4endl 00226 << " SOLID = " << sol << G4endl 00227 << " " << sol->GetName() << " of type " << sol->GetType() 00228 << G4endl; 00229 } 00230 #endif 00231 00232 //----- Check if solid exists already 00233 G4VSolid* solid = G4tgbVolumeMgr::GetInstance() 00234 ->FindG4Solid( sol->GetName() ); 00235 if( solid ) { return solid; } 00236 00237 // Give 'sol' as Boolean solids needs to call this method twice 00238 00239 #ifdef G4VERBOSE 00240 if( G4tgrMessenger::GetVerboseLevel() >= 2 ) 00241 { 00242 G4cout << " G4tgbVolume::FindOrConstructG4Solid() - " 00243 << sol->GetSolidParams().size() << G4endl; 00244 } 00245 #endif 00246 00247 std::vector<G4double> solParam; 00248 00249 // In case of BOOLEAN solids, solidParams are taken from components 00250 00251 if( sol->GetSolidParams().size() == 1) 00252 { 00253 solParam = * sol->GetSolidParams()[ 0 ]; 00254 } 00255 00256 //----------- instantiate the appropiate G4VSolid type 00257 G4String stype = sol->GetType(); 00258 G4String sname = sol->GetName(); 00259 00260 if( stype == "BOX" ) 00261 { 00262 CheckNoSolidParams( stype, 3, solParam.size() ); 00263 solid = new G4Box( sname, solParam[0], solParam[1], solParam[2] ); 00264 00265 } 00266 else if( stype == "TUBE" ) 00267 { 00268 CheckNoSolidParams( stype, 3, solParam.size() ); 00269 solid = new G4Tubs( sname, solParam[0], solParam[1], solParam[2], 00270 0.*deg, 360.*deg ); 00271 } 00272 else if( stype == "TUBS" ) 00273 { 00274 CheckNoSolidParams( stype, 5, solParam.size() ); 00275 G4double phiDelta = solParam[4]; 00276 if( std::fabs(phiDelta - twopi) < angularTolerance ) { phiDelta = twopi; } 00277 solid = new G4Tubs( sname, solParam[0], solParam[1], solParam[2], 00278 solParam[3], phiDelta ); 00279 } 00280 else if( stype == "TRAP" ) 00281 { 00282 if( solParam.size() == 11 ) 00283 { 00284 solid = new G4Trap( sname, solParam[0], solParam[1], solParam[2], 00285 solParam[3], solParam[4], solParam[5], solParam[6], 00286 solParam[7], solParam[8], solParam[9], solParam[10] ); 00287 } 00288 else if( solParam.size() == 4 ) 00289 { 00290 solid = new G4Trap( sname, solParam[0], solParam[1]/deg, 00291 solParam[2]/deg, solParam[3]); 00292 } 00293 else 00294 { 00295 G4String ErrMessage1 = "Solid type " + stype; 00296 G4String ErrMessage2 = " should have 11 or 4 parameters,\n"; 00297 G4String ErrMessage3 = "and it has " 00298 + G4UIcommand::ConvertToString(G4int(solParam.size())); 00299 G4String ErrMessage = ErrMessage1 + ErrMessage2 + ErrMessage3 + " !"; 00300 G4Exception("G4tgbVolume::FindOrConstructG4Solid()", 00301 "InvalidSetup", FatalException, ErrMessage); 00302 return 0; 00303 } 00304 00305 } 00306 else if( stype == "TRD" ) 00307 { 00308 CheckNoSolidParams( stype, 5, solParam.size() ); 00309 solid = new G4Trd( sname, solParam[0], solParam[1], solParam[2], 00310 solParam[3], solParam[4] ); 00311 } 00312 else if( stype == "PARA" ) 00313 { 00314 CheckNoSolidParams( stype, 6, solParam.size() ); 00315 solid = new G4Para( sname, solParam[0], solParam[1], solParam[2], 00316 solParam[3], solParam[4], solParam[5] ); 00317 } 00318 else if( stype == "CONE" ) 00319 { 00320 CheckNoSolidParams( stype, 5, solParam.size() ); 00321 solid = new G4Cons( sname, solParam[0], solParam[1], solParam[2], 00322 solParam[3], solParam[4], 0., 360.*deg); 00323 } 00324 else if( stype == "CONS" ) 00325 { 00326 CheckNoSolidParams( stype, 7, solParam.size() ); 00327 G4double phiDelta = solParam[6]; 00328 if( std::fabs(phiDelta - twopi) < angularTolerance ) { phiDelta = twopi; } 00329 solid = new G4Cons( sname, solParam[0], solParam[1], solParam[2], 00330 solParam[3], solParam[4], solParam[5], phiDelta); 00331 } 00332 else if( stype == "SPHERE" ) 00333 { 00334 CheckNoSolidParams( stype, 6, solParam.size() ); 00335 G4double phiDelta = solParam[3]; 00336 if( std::fabs(phiDelta - twopi) < angularTolerance ) { phiDelta = twopi; } 00337 G4double thetaDelta = solParam[5]; 00338 if( std::fabs(thetaDelta - pi) < angularTolerance ) { thetaDelta = pi; } 00339 solid = new G4Sphere( sname, solParam[0], solParam[1], solParam[2], 00340 phiDelta, solParam[4], thetaDelta); 00341 } 00342 else if( stype == "ORB" ) 00343 { 00344 CheckNoSolidParams( stype, 1, solParam.size() ); 00345 solid = new G4Orb( sname, solParam[0] ); 00346 } 00347 else if( stype == "TORUS" ) 00348 { 00349 CheckNoSolidParams( stype, 5, solParam.size() ); 00350 G4double phiDelta = solParam[4]; 00351 if( std::fabs(phiDelta - twopi) < angularTolerance ) { phiDelta = twopi; } 00352 solid = new G4Torus( sname, solParam[0], solParam[1], solParam[2], 00353 solParam[3], phiDelta ); 00354 } 00355 else if( stype == "POLYCONE" ) 00356 { 00357 size_t nplanes = size_t(solParam[2]); 00358 G4bool genericPoly = false; 00359 if( solParam.size() == 3+nplanes*3 ) 00360 { 00361 genericPoly = true; 00362 } 00363 else if( solParam.size() == 3+nplanes*2 ) 00364 { 00365 genericPoly = false; 00366 } 00367 else 00368 { 00369 G4String Err1 = "Solid type " + stype + " should have "; 00370 G4String Err2 = G4UIcommand::ConvertToString(G4int(3+nplanes*3)) 00371 + " (Z,Rmin,Rmax)\n"; 00372 G4String Err3 = "or " + G4UIcommand::ConvertToString(G4int(3+nplanes*2)); 00373 G4String Err4 = " (RZ corners) parameters,\n"; 00374 G4String Err5 = "and it has " 00375 + G4UIcommand::ConvertToString(G4int(solParam.size())); 00376 G4String ErrMessage = Err1 + Err2 + Err3 + Err4 + Err5 + " !"; 00377 G4Exception("G4tgbVolume::FindOrConstructG4Solid()", 00378 "InvalidSetup", FatalException, ErrMessage); 00379 return 0; 00380 } 00381 00382 if( genericPoly ) 00383 { 00384 std::vector<G4double>* z_p = new std::vector<G4double>; 00385 std::vector<G4double>* rmin_p = new std::vector<G4double>; 00386 std::vector<G4double>* rmax_p = new std::vector<G4double>; 00387 for( size_t ii = 0; ii < nplanes; ii++ ) 00388 { 00389 (*z_p).push_back( solParam[3+3*ii] ); 00390 (*rmin_p).push_back( solParam[3+3*ii+1] ); 00391 (*rmax_p).push_back( solParam[3+3*ii+2] ); 00392 } 00393 G4double phiTotal = solParam[1]; 00394 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; } 00395 solid = new G4Polycone( sname, solParam[0], phiTotal, // start,delta-phi 00396 nplanes, // sections 00397 &((*z_p)[0]), &((*rmin_p)[0]), &((*rmax_p)[0])); 00398 } 00399 else 00400 { 00401 std::vector<G4double>* R_c = new std::vector<G4double>; 00402 std::vector<G4double>* Z_c = new std::vector<G4double>; 00403 for( size_t ii = 0; ii < nplanes; ii++ ) 00404 { 00405 (*R_c).push_back( solParam[3+2*ii] ); 00406 (*Z_c).push_back( solParam[3+2*ii+1] ); 00407 } 00408 G4double phiTotal = solParam[1]; 00409 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; } 00410 solid = new G4Polycone( sname, solParam[0], phiTotal, // start,delta-phi 00411 nplanes, // sections 00412 &((*R_c)[0]), &((*Z_c)[0])); 00413 } 00414 00415 } 00416 else if( stype == "POLYHEDRA" ) 00417 { 00418 size_t nplanes = size_t(solParam[3]); 00419 G4bool genericPoly = false; 00420 if( solParam.size() == 4+nplanes*3 ) 00421 { 00422 genericPoly = true; 00423 } 00424 else if( solParam.size() == 4+nplanes*2 ) 00425 { 00426 genericPoly = false; 00427 } 00428 else 00429 { 00430 G4String Err1 = "Solid type " + stype + " should have "; 00431 G4String Err2 = G4UIcommand::ConvertToString(G4int(4+nplanes*3)) 00432 + " (Z,Rmin,Rmax)\n"; 00433 G4String Err3 = "or " + G4UIcommand::ConvertToString(G4int(4+nplanes*2)); 00434 G4String Err4 = " (RZ corners) parameters,\n"; 00435 G4String Err5 = "and it has " 00436 + G4UIcommand::ConvertToString(G4int(solParam.size())); 00437 G4String ErrMessage = Err1 + Err2 + Err3 + Err4 + Err5 + " !"; 00438 G4Exception("G4tgbVolume::FindOrConstructG4Solid()", 00439 "InvalidSetup", FatalException, ErrMessage); 00440 return 0; 00441 } 00442 00443 if( genericPoly ) 00444 { 00445 std::vector<G4double>* z_p = new std::vector<G4double>; 00446 std::vector<G4double>* rmin_p = new std::vector<G4double>; 00447 std::vector<G4double>* rmax_p = new std::vector<G4double>; 00448 for( size_t ii = 0; ii < nplanes; ii++ ) 00449 { 00450 (*z_p).push_back( solParam[4+3*ii] ); 00451 (*rmin_p).push_back( solParam[4+3*ii+1] ); 00452 (*rmax_p).push_back( solParam[4+3*ii+2] ); 00453 } 00454 G4double phiTotal = solParam[1]; 00455 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; } 00456 solid = new G4Polyhedra( sname, solParam[0], phiTotal, 00457 G4int(solParam[2]), nplanes, 00458 &((*z_p)[0]), &((*rmin_p)[0]), &((*rmax_p)[0])); 00459 } 00460 else 00461 { 00462 std::vector<G4double>* R_c = new std::vector<G4double>; 00463 std::vector<G4double>* Z_c = new std::vector<G4double>; 00464 for( size_t ii = 0; ii < nplanes; ii++ ) 00465 { 00466 (*R_c).push_back( solParam[4+2*ii] ); 00467 (*Z_c).push_back( solParam[4+2*ii+1] ); 00468 } 00469 G4double phiTotal = solParam[1]; 00470 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; } 00471 solid = new G4Polyhedra( sname, solParam[0], phiTotal, 00472 G4int(solParam[2]), nplanes, 00473 &((*R_c)[0]), &((*Z_c)[0])); 00474 } 00475 } 00476 else if( stype == "ELLIPTICALTUBE" ) 00477 { 00478 CheckNoSolidParams( stype, 3, solParam.size() ); 00479 solid = new G4EllipticalTube( sname, solParam[0], solParam[1], solParam[2]); 00480 } 00481 else if( stype == "ELLIPSOID" ) 00482 { 00483 CheckNoSolidParams( stype, 5, solParam.size() ); 00484 solid = new G4Ellipsoid( sname, solParam[0], solParam[1], solParam[2], 00485 solParam[3], solParam[4] ); 00486 } 00487 else if( stype == "ELLIPTICALCONE" ) 00488 { 00489 CheckNoSolidParams( stype, 4, solParam.size() ); 00490 solid = new G4EllipticalCone( sname, solParam[0], solParam[1], 00491 solParam[2], solParam[3] ); 00492 } 00493 else if( stype == "HYPE" ) 00494 { 00495 CheckNoSolidParams( stype, 5, solParam.size() ); 00496 solid = new G4Hype( sname, solParam[0], solParam[1], solParam[2], 00497 solParam[3], solParam[4] ); 00498 } 00499 else if( stype == "TET" ) 00500 { 00501 CheckNoSolidParams( stype, 12, solParam.size() ); 00502 G4ThreeVector anchor(solParam[0], solParam[1], solParam[2]); 00503 G4ThreeVector p2(solParam[3], solParam[4], solParam[5]); 00504 G4ThreeVector p3(solParam[6], solParam[7], solParam[8]); 00505 G4ThreeVector p4(solParam[9], solParam[10], solParam[11]); 00506 solid = new G4Tet( sname, anchor, p2, p3, p4 ); 00507 } 00508 else if( stype == "TWISTEDBOX" ) 00509 { 00510 CheckNoSolidParams( stype, 4, solParam.size() ); 00511 solid = new G4TwistedBox( sname, solParam[0], solParam[1], 00512 solParam[2], solParam[3]); 00513 } 00514 else if( stype == "TWISTEDTRAP" ) 00515 { 00516 CheckNoSolidParams( stype, 11, solParam.size() ); 00517 solid = new G4TwistedTrap( sname, solParam[0], solParam[1], solParam[2], 00518 solParam[3], solParam[4], solParam[5], solParam[6], 00519 solParam[7], solParam[8], solParam[9], solParam[10] ); 00520 } 00521 else if( stype == "TWISTEDTRD" ) 00522 { 00523 CheckNoSolidParams( stype, 6, solParam.size() ); 00524 solid = new G4TwistedTrd( sname, solParam[0], solParam[1], solParam[2], 00525 solParam[3], solParam[4], solParam[5]); 00526 } 00527 else if( stype == "TWISTEDTUBS" ) 00528 { 00529 CheckNoSolidParams( stype, 5, solParam.size() ); 00530 G4double phiTotal = solParam[4]; 00531 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; } 00532 solid = new G4TwistedTubs( sname, solParam[0], solParam[1], solParam[2], 00533 solParam[3], phiTotal); 00534 } 00535 else if( stype == "BREPBOX" ) // EntityType is = "Closed_Shell" 00536 { 00537 CheckNoSolidParams( stype, 24, solParam.size() ); 00538 std::vector<G4Point3D> points; 00539 for( size_t ii = 0; ii < 8; ii++ ) 00540 { 00541 points.push_back( G4Point3D(solParam[ii*3+0], 00542 solParam[ii*3+1], 00543 solParam[ii*3+2]) ); 00544 } 00545 solid = new G4BREPSolidBox( sname, points[0], points[1], points[2], 00546 points[3], points[4], points[5], points[6], 00547 points[7] ); 00548 } 00549 else if( stype == "BREPCYLINDER" ) // EntityType is = "Closed_Shell" 00550 { 00551 CheckNoSolidParams( stype, 11, solParam.size() ); 00552 solid = new G4BREPSolidCylinder( sname, 00553 G4ThreeVector( solParam[0], solParam[1], solParam[2] ), 00554 G4ThreeVector( solParam[3], solParam[4], solParam[5] ), 00555 G4ThreeVector( solParam[6], solParam[7], solParam[8] ), 00556 solParam[9], solParam[10] ); 00557 } 00558 else if( stype == "BREPCONE" ) // EntityType is = "Closed_Shell" 00559 { 00560 CheckNoSolidParams( stype, 12, solParam.size() ); 00561 solid = new G4BREPSolidCone( sname, 00562 G4ThreeVector( solParam[0], solParam[1], solParam[2] ), 00563 G4ThreeVector( solParam[3], solParam[4], solParam[5] ), 00564 G4ThreeVector( solParam[6], solParam[7], solParam[8] ), 00565 solParam[9], solParam[10], solParam[11] ); 00566 } 00567 else if( stype == "BREPSPHERE" ) // EntityType is = "Closed_Shell" 00568 { 00569 CheckNoSolidParams( stype, 10, solParam.size() ); 00570 solid = new G4BREPSolidSphere( sname, 00571 G4ThreeVector( solParam[0], solParam[1], solParam[2] ), 00572 G4ThreeVector( solParam[3], solParam[4], solParam[5] ), 00573 G4ThreeVector( solParam[6], solParam[7], solParam[8] ), 00574 solParam[9] ); 00575 00576 } 00577 else if( stype == "BREPTORUS" ) // EntityType is = "Closed_Shell" 00578 { 00579 CheckNoSolidParams( stype, 11, solParam.size() ); 00580 solid = new G4BREPSolidTorus( sname, 00581 G4ThreeVector( solParam[0], solParam[1], solParam[2] ), 00582 G4ThreeVector( solParam[3], solParam[4], solParam[5] ), 00583 G4ThreeVector( solParam[6], solParam[7], solParam[8] ), 00584 solParam[9], solParam[10] ); 00585 } 00586 else if( stype == "BREPPCONE" ) // EntityType is = "Closed_Shell" 00587 { 00588 size_t nplanes = size_t(solParam[2]); 00589 CheckNoSolidParams( stype, 4+3*nplanes, solParam.size() ); 00590 std::vector<G4double>* z_p = new std::vector<G4double>; 00591 std::vector<G4double>* rmin_p = new std::vector<G4double>; 00592 std::vector<G4double>* rmax_p = new std::vector<G4double>; 00593 for( size_t ii = 0; ii < nplanes; ii++ ) 00594 { 00595 (*z_p).push_back( solParam[4+3*ii] ); 00596 (*rmin_p).push_back( solParam[4+3*ii+1] ); 00597 (*rmax_p).push_back( solParam[4+3*ii+2] ); 00598 } 00599 G4double phiTotal = solParam[1]; 00600 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; } 00601 CheckNoSolidParams( stype, 12, solParam.size() ); 00602 solid = new G4BREPSolidPCone( sname, solParam[0], phiTotal, // start,dph 00603 nplanes, // sections 00604 solParam[3], // z_start 00605 &((*z_p)[0]), &((*rmin_p)[0]), 00606 &((*rmax_p)[0])); 00607 } 00608 else if( stype == "BREPPOLYHEDRA" ) // EntityType is = "Closed_Shell" 00609 { 00610 size_t nplanes = size_t(solParam[3]); 00611 CheckNoSolidParams( stype, 5+3*nplanes, solParam.size() ); 00612 std::vector<G4double>* z_p = new std::vector<G4double>; 00613 std::vector<G4double>* rmin_p = new std::vector<G4double>; 00614 std::vector<G4double>* rmax_p = new std::vector<G4double>; 00615 for( size_t ii = 0; ii < nplanes; ii++ ) 00616 { 00617 (*z_p).push_back( solParam[5+3*ii] ); 00618 (*rmin_p).push_back( solParam[5+3*ii+1] ); 00619 (*rmax_p).push_back( solParam[5+3*ii+2] ); 00620 } 00621 G4double phiTotal = solParam[1]; 00622 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; } 00623 CheckNoSolidParams( stype, 12, solParam.size() ); 00624 solid = new G4BREPSolidPolyhedra( sname, solParam[0], phiTotal, // start,dph 00625 G4int(solParam[2]), // sides 00626 nplanes, // sections 00627 solParam[4], // z_start 00628 &((*z_p)[0]), &((*rmin_p)[0]), 00629 &((*rmax_p)[0])); 00630 } 00631 else if( stype == "BREPOPENPCONE" ) // EntityType is = "Closed_Shell" 00632 { 00633 size_t nplanes = size_t(solParam[2]); 00634 std::vector<G4double>* z_p = new std::vector<G4double>; 00635 std::vector<G4double>* rmin_p = new std::vector<G4double>; 00636 std::vector<G4double>* rmax_p = new std::vector<G4double>; 00637 for( size_t ii = 0; ii < nplanes; ii++ ) 00638 { 00639 (*z_p).push_back( solParam[4+3*ii] ); 00640 (*rmin_p).push_back( solParam[4+3*ii+1] ); 00641 (*rmax_p).push_back( solParam[4+3*ii+2] ); 00642 } 00643 G4double phiTotal = solParam[1]; 00644 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; } 00645 CheckNoSolidParams( stype, 12, solParam.size() ); 00646 solid = new G4BREPSolidOpenPCone( sname, solParam[0], phiTotal, // start,dph 00647 nplanes, // sections 00648 solParam[3], // z_start 00649 &((*z_p)[0]), &((*rmin_p)[0]), 00650 &((*rmax_p)[0])); 00651 } 00652 else if( stype == "TESSELLATED" ) 00653 { 00654 G4int nFacets = G4int(solParam[0]); 00655 G4int jj = 0; 00656 solid = new G4TessellatedSolid(sname); 00657 G4TessellatedSolid* solidTS = (G4TessellatedSolid*)(solid); 00658 G4VFacet* facet=0; 00659 00660 for( G4int ii = 0; ii < nFacets; ii++){ 00661 G4int nPoints = G4int(solParam[jj+1]); 00662 if( G4int(solParam.size()) < jj + nPoints*3 + 2 ) { 00663 G4String Err1 = "Too small number of parameters in tesselated solid, it should be at least " + G4UIcommand::ConvertToString(jj + nPoints*3 + 2 ); 00664 G4String Err2 = " facet number " + G4UIcommand::ConvertToString(ii); 00665 G4String Err3 = " number of parameters is " + G4UIcommand::ConvertToString(G4int(solParam.size())); 00666 G4String ErrMessage = Err1 + Err2 + Err3 + " !"; 00667 G4Exception("G4tgbVolume::FindOrConstructG4Solid()", 00668 "InvalidSetup", FatalException, ErrMessage); 00669 return 0; 00670 } 00671 00672 if( nPoints == 3 ) 00673 { 00674 G4ThreeVector pt0(solParam[jj+2],solParam[jj+3],solParam[jj+4]); 00675 G4ThreeVector vt1(solParam[jj+5],solParam[jj+6],solParam[jj+7]); 00676 G4ThreeVector vt2(solParam[jj+8],solParam[jj+9],solParam[jj+10]); 00677 G4FacetVertexType vertexType = ABSOLUTE; 00678 if( solParam[jj+11] == 0 ) 00679 { 00680 vertexType = ABSOLUTE; 00681 } 00682 else if( solParam[jj+11] == 1 ) 00683 { 00684 vertexType = RELATIVE; 00685 } 00686 else 00687 { 00688 G4String Err1 = "Wrong number of vertex type in tesselated solid, it should be 0 =ABSOLUTE) or 1 (=RELATIVE)"; 00689 G4String Err2 = " facet number " + G4UIcommand::ConvertToString(G4int(ii)); 00690 G4String Err3 = " vertex type is " + G4UIcommand::ConvertToString(solParam[jj+11]); 00691 G4String ErrMessage = Err1 + Err2 + Err3 + " !"; 00692 G4Exception("G4tgbVolume::FindOrConstructG4Solid()", 00693 "InvalidSetup", FatalException, ErrMessage); 00694 return 0; 00695 } 00696 facet = new G4TriangularFacet( pt0, vt1, vt2, vertexType ); 00697 } 00698 else if( nPoints == 4 ) 00699 { 00700 G4ThreeVector pt0(solParam[jj+2],solParam[jj+3],solParam[jj+4]); 00701 G4ThreeVector vt1(solParam[jj+5],solParam[jj+6],solParam[jj+7]); 00702 G4ThreeVector vt2(solParam[jj+8],solParam[jj+9],solParam[jj+10]); 00703 G4ThreeVector vt3(solParam[jj+11],solParam[jj+12],solParam[jj+13]); 00704 G4FacetVertexType vertexType = ABSOLUTE; 00705 if( solParam[jj+14] == 0 ) 00706 { 00707 vertexType = ABSOLUTE; 00708 } 00709 else if( solParam[jj+14] == 1 ) 00710 { 00711 vertexType = RELATIVE; 00712 } 00713 else 00714 { 00715 G4String Err1 = "Wrong number of vertex type in tesselated solid, it should be 0 =ABSOLUTE) or 1 (=RELATIVE)"; 00716 G4String Err2 = " facet number " + G4UIcommand::ConvertToString(G4int(ii)); 00717 G4String Err3 = " vertex type is " + G4UIcommand::ConvertToString(solParam[jj+14]); 00718 G4String ErrMessage = Err1 + Err2 + Err3 + " !"; 00719 G4Exception("G4tgbVolume::FindOrConstructG4Solid()", 00720 "InvalidSetup", FatalException, ErrMessage); 00721 return 0; 00722 } 00723 facet = new G4QuadrangularFacet( pt0, vt1, vt2, vt3, vertexType ); 00724 } 00725 else 00726 { 00727 G4String Err1 = "Wrong number of points in tesselated solid, it should be 3 or 4"; 00728 G4String Err2 = " facet number " + G4UIcommand::ConvertToString(G4int(ii)); 00729 G4String Err3 = " number of points is " + G4UIcommand::ConvertToString(G4int(nPoints)); 00730 G4String ErrMessage = Err1 + Err2 + Err3 + " !"; 00731 G4Exception("G4tgbVolume::FindOrConstructG4Solid()", 00732 "InvalidSetup", FatalException, ErrMessage); 00733 return 0; 00734 } 00735 00736 solidTS->AddFacet( facet ); 00737 jj += nPoints*3 + 2; 00738 } 00739 00740 } 00741 else if( stype == "EXTRUDED" ) 00742 { 00743 std::vector<G4TwoVector> polygonList; 00744 std::vector<G4ExtrudedSolid::ZSection> zsectionList; 00745 G4int nPolygons = G4int(solParam[0]); 00746 G4int ii = 1; 00747 G4int nMax = nPolygons*2+1; 00748 for( ;ii < nMax; ii+=2 ) 00749 { 00750 polygonList.push_back( G4TwoVector(solParam[ii],solParam[ii+1]) ); 00751 } 00752 G4int nZSections = G4int(solParam[ii]); 00753 nMax = nPolygons*2 + nZSections*4 + 2; 00754 ii++; 00755 for( ; ii < nMax; ii+=4 ) 00756 { 00757 G4TwoVector offset(solParam[ii+1],solParam[ii+2]); 00758 zsectionList.push_back( G4ExtrudedSolid::ZSection(solParam[ii],offset,solParam[ii+3]) ); 00759 } 00760 solid = new G4ExtrudedSolid( sname, polygonList, zsectionList ); 00761 00762 } 00763 else if( stype.substr(0,7) == "Boolean" ) 00764 { 00765 const G4tgrSolidBoolean* solb = dynamic_cast<const G4tgrSolidBoolean*>(sol); 00766 if (!solb) 00767 { 00768 G4Exception("G4tgbVolume::FindOrConstructG4Solid()", 00769 "InvalidSetup", FatalException, "Invalid Solid pointer"); 00770 return 0; 00771 } 00772 G4VSolid* sol1 = FindOrConstructG4Solid( solb->GetSolid(0)); 00773 G4VSolid* sol2 = FindOrConstructG4Solid( solb->GetSolid(1)); 00774 G4RotationMatrix* relRotMat = G4tgbRotationMatrixMgr::GetInstance() 00775 ->FindOrBuildG4RotMatrix( sol->GetRelativeRotMatName() ); 00776 G4ThreeVector relPlace = solb->GetRelativePlace(); 00777 00778 if( stype == "Boolean_UNION" ) 00779 { 00780 solid = new G4UnionSolid( sname, sol1, sol2, relRotMat, relPlace ); 00781 } 00782 else if( stype == "Boolean_SUBTRACTION" ) 00783 { 00784 solid = new G4SubtractionSolid( sname, sol1, sol2, relRotMat, relPlace ); 00785 } 00786 else if( stype == "Boolean_INTERSECTION" ) 00787 { 00788 solid = new G4IntersectionSolid( sname, sol1, sol2, relRotMat, relPlace ); 00789 } 00790 else 00791 { 00792 G4String ErrMessage = "Unknown Boolean type " + stype; 00793 G4Exception("G4tgbVolume::FindOrConstructG4Solid()", 00794 "InvalidSetup", FatalException, ErrMessage); 00795 return 0; 00796 } 00797 } 00798 else 00799 { 00800 G4String ErrMessage = "Solids of type " + stype 00801 + " not implemented yet, sorry..."; 00802 G4Exception("G4tgbVolume::FindOrConstructG4Solid()", "NotImplemented", 00803 FatalException, ErrMessage); 00804 return 0; 00805 } 00806 00807 #ifdef G4VERBOSE 00808 if( G4tgrMessenger::GetVerboseLevel() >= 2 ) 00809 { 00810 G4cout << " G4tgbVolume::FindOrConstructG4Solid()" << G4endl 00811 << " Created solid " << sname 00812 << " of type " << solid->GetEntityType() << G4endl; 00813 } 00814 #endif 00815 00816 #ifdef G4VERBOSE 00817 if( G4tgrMessenger::GetVerboseLevel() >= 1 ) 00818 { 00819 G4cout << " Constructing new G4Solid: " 00820 << *solid << G4endl; 00821 } 00822 #endif 00823 00824 return solid; 00825 }
const G4double* G4tgbVolume::GetColour | ( | ) | const [inline] |
Definition at line 110 of file G4tgbVolume.hh.
References G4tgrVolume::GetColour().
Referenced by ConstructG4LogVol().
00110 { return theTgrVolume->GetColour(); }
const G4String& G4tgbVolume::GetName | ( | ) | const [inline] |
Definition at line 108 of file G4tgbVolume.hh.
References G4tgrVolume::GetName().
Referenced by BuildSolidForDivision(), ConstructG4LogVol(), ConstructG4PhysVol(), ConstructG4Volumes(), and G4tgbVolumeMgr::RegisterMe().
00108 { return theTgrVolume->GetName(); }
G4bool G4tgbVolume::GetVisibility | ( | ) | const [inline] |
Definition at line 109 of file G4tgbVolume.hh.
References G4tgrVolume::GetVisibility().
Referenced by ConstructG4LogVol().
00109 { return theTgrVolume->GetVisibility(); }
void G4tgbVolume::SetCutsInEnergy | ( | G4LogicalVolume * | logvol, | |
std::map< G4String, G4double > | cuts | |||
) |
void G4tgbVolume::SetCutsInRange | ( | G4LogicalVolume * | logvol, | |
std::map< G4String, G4double > | cuts | |||
) |