Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4tgbVolume Class Reference

#include <G4tgbVolume.hh>

Public Member Functions

 G4tgbVolume ()
 
 ~G4tgbVolume ()
 
 G4tgbVolume (G4tgrVolume *vol)
 
void ConstructG4Volumes (const G4tgrPlace *place, const G4LogicalVolume *parentLV)
 
G4VSolidFindOrConstructG4Solid (const G4tgrSolid *vol)
 
G4LogicalVolumeConstructG4LogVol (const G4VSolid *solid)
 
G4VPhysicalVolumeConstructG4PhysVol (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)
 
G4VSolidBuildSolidForDivision (G4VSolid *parentSolid, EAxis axis)
 
const G4StringGetName () const
 
G4bool GetVisibility () const
 
const G4doubleGetColour () const
 

Detailed Description

Definition at line 66 of file G4tgbVolume.hh.

Constructor & Destructor Documentation

◆ G4tgbVolume() [1/2]

G4tgbVolume::G4tgbVolume ( )

Definition at line 116 of file G4tgbVolume.cc.

117 : theTgrVolume(0), theG4AssemblyVolume(0)
118{
119}

◆ ~G4tgbVolume()

G4tgbVolume::~G4tgbVolume ( )

Definition at line 123 of file G4tgbVolume.cc.

124{
125}

◆ G4tgbVolume() [2/2]

G4tgbVolume::G4tgbVolume ( G4tgrVolume vol)

Definition at line 129 of file G4tgbVolume.cc.

130{
131 theTgrVolume = vol;
132 theG4AssemblyVolume = 0;
133}

Member Function Documentation

◆ BuildSolidForDivision()

G4VSolid * G4tgbVolume::BuildSolidForDivision ( G4VSolid parentSolid,
EAxis  axis 
)

Definition at line 1310 of file G4tgbVolume.cc.

1311{
1312 G4VSolid* solid=0;
1313 G4double redf = (parentSolid->GetExtent().GetXmax()-parentSolid->GetExtent().GetXmin());
1314 redf = std::min(redf,parentSolid->GetExtent().GetYmax()-parentSolid->GetExtent().GetYmin());
1315 redf = std::min(redf,parentSolid->GetExtent().GetZmax()-parentSolid->GetExtent().GetZmin());
1316 redf *= 0.001; //make daugther much smaller, to fit in parent
1317
1318 if( parentSolid->GetEntityType() == "G4Box" )
1319 {
1320 G4Box* psolid = (G4Box*)(parentSolid);
1321 solid = new G4Box(GetName(), psolid->GetXHalfLength()*redf,
1322 psolid->GetZHalfLength()*redf,
1323 psolid->GetZHalfLength()*redf);
1324 }
1325 else if ( parentSolid->GetEntityType() == "G4Tubs" )
1326 {
1327 G4Tubs* psolid = (G4Tubs*)(parentSolid);
1328 solid = new G4Tubs( GetName(), psolid->GetInnerRadius()*redf,
1329 psolid->GetOuterRadius()*redf,
1330 psolid->GetZHalfLength()*redf,
1331 psolid->GetSPhi(), psolid->GetDPhi());
1332 }
1333 else if ( parentSolid->GetEntityType() == "G4Cons" )
1334 {
1335 G4Cons* psolid = (G4Cons*)(parentSolid);
1336 solid = new G4Cons( GetName(), psolid->GetInnerRadiusMinusZ()*redf,
1337 psolid->GetOuterRadiusMinusZ()*redf,
1338 psolid->GetInnerRadiusPlusZ()*redf,
1339 psolid->GetOuterRadiusPlusZ()*redf,
1340 psolid->GetZHalfLength()*redf,
1341 psolid->GetSPhi(), psolid->GetDPhi());
1342 }
1343 else if ( parentSolid->GetEntityType() == "G4Trd" )
1344 {
1345 G4Trd* psolid = (G4Trd*)(parentSolid);
1346 G4double mpDx1 = psolid->GetXHalfLength1();
1347 G4double mpDx2 = psolid->GetXHalfLength2();
1348
1349 if( axis == kXAxis && std::fabs(mpDx1 - mpDx2) > G4GeometryTolerance::GetInstance()->GetSurfaceTolerance() )
1350 {
1351 solid = new G4Trap( GetName(), psolid->GetZHalfLength()*redf,
1352 psolid->GetYHalfLength1()*redf,
1353 psolid->GetXHalfLength2()*redf,
1354 psolid->GetXHalfLength1()*redf );
1355 }
1356 else
1357 {
1358 solid = new G4Trd( GetName(), psolid->GetXHalfLength1()*redf,
1359 psolid->GetXHalfLength2()*redf,
1360 psolid->GetYHalfLength1()*redf,
1361 psolid->GetYHalfLength2()*redf,
1362 psolid->GetZHalfLength()*redf);
1363 }
1364
1365 }
1366 else if ( parentSolid->GetEntityType() == "G4Para" )
1367 {
1368 G4Para* psolid = (G4Para*)(parentSolid);
1369 solid = new G4Para( GetName(), psolid->GetXHalfLength()*redf,
1370 psolid->GetYHalfLength()*redf,
1371 psolid->GetZHalfLength()*redf,
1372 std::atan(psolid->GetTanAlpha()),
1373 psolid->GetSymAxis().theta(),
1374 psolid->GetSymAxis().phi() );
1375 }
1376 else if ( parentSolid->GetEntityType() == "G4Polycone" )
1377 {
1378 G4Polycone* psolid = (G4Polycone*)(parentSolid);
1379 G4PolyconeHistorical origParam = *(psolid->GetOriginalParameters());
1380 for( G4int ii = 0; ii < origParam.Num_z_planes; ii++ )
1381 {
1382 origParam.Rmin[ii] = origParam.Rmin[ii]*redf;
1383 origParam.Rmax[ii] = origParam.Rmax[ii]*redf;
1384 }
1385 solid = new G4Polycone( GetName(), psolid->GetStartPhi(),
1386 psolid->GetEndPhi(),
1387 origParam.Num_z_planes, origParam.Z_values,
1388 origParam.Rmin, origParam.Rmax);
1389
1390 }
1391 else if ( parentSolid->GetEntityType() == "G4Polyhedra" )
1392 {
1393 G4Polyhedra* psolid = (G4Polyhedra*)(parentSolid);
1394 G4PolyhedraHistorical origParam = *(psolid->GetOriginalParameters());
1395 for( G4int ii = 0; ii < origParam.Num_z_planes; ii++ )
1396 {
1397 origParam.Rmin[ii] = origParam.Rmin[ii]*redf;
1398 origParam.Rmax[ii] = origParam.Rmax[ii]*redf;
1399 }
1400 solid = new G4Polyhedra( GetName(), psolid->GetStartPhi(),
1401 psolid->GetEndPhi(),
1402 psolid->GetNumSide(),
1403 origParam.Num_z_planes, origParam.Z_values,
1404 origParam.Rmin, origParam.Rmax);
1405 }
1406 else
1407 {
1408 G4String ErrMessage = "Solid type not supported. VOLUME= " + GetName()
1409 + " Solid type= " + parentSolid->GetEntityType() + "\n"
1410 + "Only supported types are: G4Box, G4Tubs, G4Cons,"
1411 + " G4Trd, G4Para, G4Polycone, G4Polyhedra.";
1412 G4Exception("G4tgbVolume::BuildSolidForDivision()", "NotImplemented",
1413 FatalException, ErrMessage);
1414 return 0;
1415 }
1416
1417#ifdef G4VERBOSE
1419 {
1420 G4cout << " Constructing new G4Solid for division: "
1421 << *solid << G4endl;
1422 }
1423#endif
1424 return solid;
1425}
@ FatalException
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
double phi() const
double theta() const
Definition: G4Box.hh:55
G4double GetZHalfLength() const
G4double GetXHalfLength() const
Definition: G4Cons.hh:75
G4double GetOuterRadiusPlusZ() const
G4double GetInnerRadiusMinusZ() const
G4double GetInnerRadiusPlusZ() const
G4double GetSPhi() const
G4double GetOuterRadiusMinusZ() const
G4double GetZHalfLength() const
G4double GetDPhi() const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
Definition: G4Para.hh:77
G4double GetTanAlpha() const
G4ThreeVector GetSymAxis() const
G4double GetYHalfLength() const
G4double GetZHalfLength() const
G4double GetXHalfLength() const
G4double * Z_values
Definition: G4Polycone.hh:80
G4double GetEndPhi() const
G4double GetStartPhi() const
G4PolyconeHistorical * GetOriginalParameters() const
G4double GetEndPhi() const
G4int GetNumSide() const
G4PolyhedraHistorical * GetOriginalParameters() const
G4double GetStartPhi() const
Definition: G4Trd.hh:63
G4double GetXHalfLength2() const
G4double GetYHalfLength2() const
G4double GetXHalfLength1() const
G4double GetYHalfLength1() const
G4double GetZHalfLength() const
Definition: G4Tubs.hh:77
G4double GetZHalfLength() const
G4double GetSPhi() const
G4double GetInnerRadius() const
G4double GetOuterRadius() const
G4double GetDPhi() const
virtual G4VisExtent GetExtent() const
Definition: G4VSolid.cc:619
virtual G4GeometryType GetEntityType() const =0
G4double GetYmin() const
Definition: G4VisExtent.hh:91
G4double GetXmax() const
Definition: G4VisExtent.hh:90
G4double GetYmax() const
Definition: G4VisExtent.hh:92
G4double GetZmax() const
Definition: G4VisExtent.hh:94
G4double GetZmin() const
Definition: G4VisExtent.hh:93
G4double GetXmin() const
Definition: G4VisExtent.hh:89
const G4String & GetName() const
Definition: G4tgbVolume.hh:108
static G4int GetVerboseLevel()
@ kXAxis
Definition: geomdefs.hh:54
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Referenced by ConstructG4PhysVol().

◆ CheckNoSolidParams()

void G4tgbVolume::CheckNoSolidParams ( const G4String solidType,
const unsigned int  NoParamExpected,
const unsigned int  NoParam 
)

Definition at line 832 of file G4tgbVolume.cc.

835{
836 if( NoParamExpected != NoParam )
837 {
838 G4String Err1 = "Solid type " + solidType + " should have ";
839 G4String Err2 = G4UIcommand::ConvertToString(G4int(NoParamExpected))
840 + " parameters,\n";
841 G4String Err3 = "and it has "
843 G4String ErrMessage = Err1 + Err2 + Err3 + " !";
844 G4Exception("G4tgbVolume::CheckNoSolidParams()", "InvalidSetup",
845 FatalException, ErrMessage);
846 }
847}
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:349

Referenced by FindOrConstructG4Solid().

◆ ConstructG4LogVol()

G4LogicalVolume * G4tgbVolume::ConstructG4LogVol ( const G4VSolid solid)

Definition at line 851 of file G4tgbVolume.cc.

852{
853 G4LogicalVolume* logvol;
854
855#ifdef G4VERBOSE
857 {
858 G4cout << " G4tgbVolume::ConstructG4LogVol() - " << GetName() << G4endl;
859 }
860#endif
861
862 //----------- Get the material first
864 ->FindOrBuildG4Material( theTgrVolume->GetMaterialName() );
865 if( mate == 0 )
866 {
867 G4String ErrMessage = "Material not found "
868 + theTgrVolume->GetMaterialName()
869 + " for volume " + GetName() + ".";
870 G4Exception("G4tgbVolume::ConstructG4LogVol()", "InvalidSetup",
871 FatalException, ErrMessage);
872 }
873#ifdef G4VERBOSE
875 {
876 G4cout << " G4tgbVolume::ConstructG4LogVol() -"
877 << " Material constructed: " << mate->GetName() << G4endl;
878 }
879#endif
880
881 //---------- Construct the LV
882 logvol = new G4LogicalVolume( const_cast<G4VSolid*>(solid),
883 const_cast<G4Material*>(mate), GetName() );
884
885#ifdef G4VERBOSE
887 {
888 G4cout << " Constructing new G4LogicalVolume: "
889 << logvol->GetName() << " mate " << mate->GetName() << G4endl;
890 }
891#endif
892
893 //---------- Set visibility and colour
894 if( !GetVisibility() || GetColour()[0] != -1 )
895 {
896 G4VisAttributes* visAtt = new G4VisAttributes();
897#ifdef G4VERBOSE
899 {
900 G4cout << " Constructing new G4VisAttributes: "
901 << *visAtt << G4endl;
902 }
903#endif
904
905 if( !GetVisibility() )
906 {
907 visAtt->SetVisibility( GetVisibility() );
908 }
909 else if( GetColour()[0] != -1 )
910 {
911 // this else should not be necessary, because if the visibility
912 // is set to off, colour should have no effect. But it does not
913 // work: if you set colour and vis off, it is visualized!?!?!?
914
915 const G4double* col = GetColour();
916 if( col[3] == -1. )
917 {
918 visAtt->SetColour( G4Colour(col[0],col[1],col[2]));
919 }
920 else
921 {
922 visAtt->SetColour( G4Colour(col[0],col[1],col[2],col[3]));
923 }
924 }
925 logvol->SetVisAttributes(visAtt);
926 }
927
928#ifdef G4VERBOSE
930 {
931 G4cout << " G4tgbVolume::ConstructG4LogVol() -"
932 << " Created logical volume: " << GetName() << G4endl;
933 }
934#endif
935
936 return logvol;
937}
G4String GetName() const
void SetVisAttributes(const G4VisAttributes *pVA)
const G4String & GetName() const
Definition: G4Material.hh:177
void SetVisibility(G4bool)
void SetColour(const G4Colour &)
G4Material * FindOrBuildG4Material(const G4String &name, G4bool bMustExist=1)
static G4tgbMaterialMgr * GetInstance()
const G4double * GetColour() const
Definition: G4tgbVolume.hh:110
G4bool GetVisibility() const
Definition: G4tgbVolume.hh:109
const G4String & GetMaterialName() const
Definition: G4tgrVolume.hh:93

Referenced by ConstructG4Volumes().

◆ ConstructG4PhysVol()

G4VPhysicalVolume * G4tgbVolume::ConstructG4PhysVol ( const G4tgrPlace place,
const G4LogicalVolume currentLV,
const G4LogicalVolume parentLV 
)

Definition at line 942 of file G4tgbVolume.cc.

945{
946 G4VPhysicalVolume* physvol = 0;
947 G4int copyNo;
948
949 //----- Case of placement of top volume
950 if( place == 0 )
951 {
952#ifdef G4VERBOSE
954 {
955 G4cout << " G4tgbVolume::ConstructG4PhysVol() - World: "
956 << GetName() << G4endl;
957 }
958#endif
959 physvol = new G4PVPlacement(0, G4ThreeVector(),
960 const_cast<G4LogicalVolume*>(currentLV),
961 GetName(), 0, false, 0,
962 theTgrVolume->GetCheckOverlaps());
963#ifdef G4VERBOSE
965 {
966 G4cout << " Constructing new : G4PVPlacement "
967 << physvol->GetName() << G4endl;
968 }
969#endif
970 }
971 else
972 {
973 copyNo = place->GetCopyNo();
974
975#ifdef G4VERBOSE
977 {
978 G4cout << " G4tgbVolume::ConstructG4PhysVol() - " << GetName() << G4endl
979 << " inside " << parentLV->GetName() << " copy No: " << copyNo
980 << " of type= " << theTgrVolume->GetType() << G4endl
981 << " placement type= " << place->GetType() << G4endl;
982 }
983#endif
984
985 if( theTgrVolume->GetType() == "VOLSimple" )
986 {
987 //----- Get placement
988#ifdef G4VERBOSE
990 {
991 G4cout << " G4tgbVolume::ConstructG4PhysVol() - Placement type = "
992 << place->GetType() << G4endl;
993 }
994#endif
995
996 //--------------- If it is G4tgrPlaceSimple
997 if( place->GetType() == "PlaceSimple" )
998 {
999 //----- Get rotation matrix
1000 G4tgrPlaceSimple* placeSimple = (G4tgrPlaceSimple*)place;
1001 G4String rmName = placeSimple->GetRotMatName();
1002
1004 ->FindOrBuildG4RotMatrix( rmName );
1005 //----- Place volume in mother
1006 G4double check = (rotmat->colX().cross(rotmat->colY()))*rotmat->colZ();
1007 G4double tol = 1.0e-3;
1008 //---- Check that matrix is ortogonal
1009 if (1-std::abs(check)>tol)
1010 {
1011 G4cerr << " Matrix : " << rmName << " " << rotmat->colX()
1012 << " " << rotmat->colY() << " " << rotmat->colZ() << G4endl
1013 << " product x X y * z = " << check << " x X y "
1014 << rotmat->colX().cross(rotmat->colY()) << G4endl;
1015 G4String ErrMessage = "Rotation is not ortogonal " + rmName + " !";
1016 G4Exception("G4tgbVolume::ConstructG4PhysVol()",
1017 "InvalidSetup", FatalException, ErrMessage);
1018 //---- Check if it is reflection
1019 }
1020 else if (1+check<=tol)
1021 {
1022 G4Translate3D transl = place->GetPlacement();
1023 G4Transform3D trfrm = transl * G4Rotate3D(*rotmat);
1024 physvol = (G4ReflectionFactory::Instance()->Place(trfrm, GetName(),
1025 const_cast<G4LogicalVolume*>(currentLV),
1026 const_cast<G4LogicalVolume*>(parentLV),
1027 false, copyNo, false )).first;
1028 }
1029 else
1030 {
1031#ifdef G4VERBOSE
1033 {
1034 G4cout << "Construction new G4VPhysicalVolume"
1035 << " through G4ReflectionFactory " << GetName()
1036 << " in volume " << parentLV->GetName()
1037 << " copyNo " << copyNo
1038 << " position " << place->GetPlacement()
1039 << " ROT " << rotmat->colX()
1040 << " " << rotmat->colY()
1041 << " " << rotmat->colZ() << G4endl;
1042 }
1043#endif
1044 physvol = new G4PVPlacement( rotmat, place->GetPlacement(),
1045 const_cast<G4LogicalVolume*>(currentLV),
1046 GetName(),
1047 const_cast<G4LogicalVolume*>(parentLV),
1048 false, copyNo,
1049 theTgrVolume->GetCheckOverlaps());
1050 }
1051
1052 //--------------- If it is G4tgrPlaceParam
1053 }
1054 else if( place->GetType() == "PlaceParam" )
1055 {
1057
1058 //----- See what parameterisation type
1059#ifdef G4VERBOSE
1061 {
1062 G4cout << " G4tgbVolume::ConstructG4PhysVol() -" << G4endl
1063 << " param: " << GetName() << " in " << parentLV->GetName()
1064 << " param type= " << dp->GetParamType() << G4endl;
1065 }
1066#endif
1067
1069
1070 if( (dp->GetParamType() == "CIRCLE")
1071 || (dp->GetParamType() == "CIRCLE_XY")
1072 || (dp->GetParamType() == "CIRCLE_XZ")
1073 || (dp->GetParamType() == "CIRCLE_YZ") )
1074 {
1075 param = new G4tgbPlaceParamCircle(dp);
1076
1077 }
1078 else if( (dp->GetParamType() == "LINEAR")
1079 || (dp->GetParamType() == "LINEAR_X")
1080 || (dp->GetParamType() == "LINEAR_Y")
1081 || (dp->GetParamType() == "LINEAR_Z") )
1082 {
1083 param = new G4tgbPlaceParamLinear(dp);
1084
1085 }
1086 else if( (dp->GetParamType() == "SQUARE")
1087 || (dp->GetParamType() == "SQUARE_XY")
1088 || (dp->GetParamType() == "SQUARE_XZ")
1089 || (dp->GetParamType() == "SQUARE_YZ") )
1090 {
1091 param = new G4tgbPlaceParamSquare(dp);
1092 }
1093 else
1094 {
1095 G4String ErrMessage = "Parameterisation has wrong type, TYPE: "
1096 + G4String(dp->GetParamType()) + " !";
1097 G4Exception("G4tgbVolume::ConstructG4PhysVol", "WrongArgument",
1098 FatalException, ErrMessage);
1099 return 0;
1100 }
1101#ifdef G4VERBOSE
1103 {
1104 G4cout << " G4tgbVolume::ConstructG4PhysVol() -" << G4endl
1105 << " New G4PVParameterised: " << GetName() << " vol "
1106 << currentLV->GetName() << " in vol " << parentLV->GetName()
1107 << " axis " << param->GetAxis() << " nCopies "
1108 << param->GetNCopies() << G4endl;
1109 }
1110#endif
1111 physvol = new G4PVParameterised(GetName(),
1112 const_cast<G4LogicalVolume*>(currentLV),
1113 const_cast<G4LogicalVolume*>(parentLV),
1114 EAxis(param->GetAxis()),
1115 param->GetNCopies(), param);
1116#ifdef G4VERBOSE
1118 {
1119 G4cout << " Constructing new G4PVParameterised: "
1120 << physvol->GetName() << " in volume " << parentLV->GetName()
1121 << " N copies " << param->GetNCopies()
1122 << " axis " << param->GetAxis() << G4endl;
1123 }
1124#endif
1125
1126 }
1127 else if( place->GetType() == "PlaceReplica" )
1128 {
1129 //--------------- If it is PlaceReplica
1130 G4tgrPlaceDivRep* dpr = (G4tgrPlaceDivRep*)place;
1131
1132#ifdef G4VERBOSE
1134 {
1135 G4cout << " G4tgbVolume::ConstructG4PhysVol() -" << G4endl
1136 << " replica" << " " << currentLV->GetName()
1137 << " in " << parentLV->GetName()
1138 << " NDiv " << dpr->GetNDiv() << " Width " << dpr->GetWidth()
1139 << " offset " << dpr->GetOffset() << G4endl;
1140 }
1141#endif
1142 physvol = new G4PVReplica(GetName(),
1143 const_cast<G4LogicalVolume*>(currentLV),
1144 const_cast<G4LogicalVolume*>(parentLV),
1145 EAxis(dpr->GetAxis()), dpr->GetNDiv(),
1146 dpr->GetWidth(), dpr->GetOffset());
1147#ifdef G4VERBOSE
1149 {
1150 G4cout << " Constructing new G4PVReplica: "
1151 << currentLV->GetName()
1152 << " in " << parentLV->GetName()
1153 << " NDiv " << dpr->GetNDiv() << " Width " << dpr->GetWidth()
1154 << " offset " << dpr->GetOffset() << G4endl;
1155 }
1156#endif
1157 }
1158 }
1159 else if( theTgrVolume->GetType() == "VOLDivision" )
1160 {
1161 G4tgrVolumeDivision* volr = (G4tgrVolumeDivision*)theTgrVolume;
1162 G4tgrPlaceDivRep* placeDiv = volr->GetPlaceDivision() ;
1163 G4VSolid* solid = BuildSolidForDivision( parentLV->GetSolid(), placeDiv->GetAxis() );
1165 ->FindOrBuildG4Material( theTgrVolume->GetMaterialName() );
1166 G4LogicalVolume* divLV = new G4LogicalVolume(solid,
1167 const_cast<G4Material*>(mate),
1168 GetName() );
1169#ifdef G4VERBOSE
1171 {
1172 G4cout << " Constructed new G4LogicalVolume for division: "
1173 << divLV->GetName() << " mate " << mate->GetName() << G4endl;
1174 }
1175#endif
1176
1177 G4DivType divType = placeDiv->GetDivType();
1178 switch (divType)
1179 {
1180 case DivByNdiv:
1181 physvol = new G4PVDivision(GetName(), (G4LogicalVolume*)divLV,
1182 const_cast<G4LogicalVolume*>(parentLV),
1183 placeDiv->GetAxis(), placeDiv->GetNDiv(),
1184 placeDiv->GetOffset());
1185#ifdef G4VERBOSE
1187 {
1188 G4cout << " Constructing new G4PVDivision by number of divisions: "
1189 << GetName() << " in " << parentLV->GetName()
1190 << " axis " << placeDiv->GetAxis()
1191 << " Ndiv " << placeDiv->GetNDiv()
1192 << " offset " << placeDiv->GetOffset() << G4endl;
1193 }
1194#endif
1195 break;
1196 case DivByWidth:
1197 physvol = new G4PVDivision(GetName(), (G4LogicalVolume*)divLV,
1198 const_cast<G4LogicalVolume*>(parentLV),
1199 placeDiv->GetAxis(), placeDiv->GetWidth(),
1200 placeDiv->GetOffset());
1201#ifdef G4VERBOSE
1203 {
1204 G4cout << " Constructing new G4PVDivision by width: "
1205 << GetName() << " in " << parentLV->GetName()
1206 << " axis " << placeDiv->GetAxis()
1207 << " width " << placeDiv->GetWidth()
1208 << " offset " << placeDiv->GetOffset() << G4endl;
1209 }
1210#endif
1211 break;
1212 case DivByNdivAndWidth:
1213 physvol = new G4PVDivision(GetName(), (G4LogicalVolume*)divLV,
1214 const_cast<G4LogicalVolume*>(parentLV),
1215 placeDiv->GetAxis(), placeDiv->GetNDiv(),
1216 placeDiv->GetWidth(),
1217 placeDiv->GetOffset());
1218#ifdef G4VERBOSE
1220 {
1221 G4cout << " Constructing new G4PVDivision by width"
1222 << " and number of divisions: "
1223 << GetName() << " in " << parentLV->GetName()
1224 << " axis " << placeDiv->GetAxis()
1225 << " Ndiv " << placeDiv->GetNDiv()
1226 << " width " << placeDiv->GetWidth()
1227 << " offset " << placeDiv->GetOffset() << G4endl;
1228 }
1229#endif
1230 break;
1231 }
1232 }
1233 else if( theTgrVolume->GetType() == "VOLAssembly" )
1234 {
1235 // Define one layer as one assembly volume
1236 G4tgrVolumeAssembly * tgrAssembly = (G4tgrVolumeAssembly *)theTgrVolume;
1237
1238 if( !theG4AssemblyVolume )
1239 {
1240 theG4AssemblyVolume = new G4AssemblyVolume;
1241#ifdef G4VERBOSE
1243 {
1244 G4cout << " Constructing new G4AssemblyVolume: "
1245 << " number of assembly components "
1246 << tgrAssembly->GetNoComponents() << G4endl;
1247 }
1248#endif
1250 for( G4int ii = 0; ii < tgrAssembly->GetNoComponents(); ii++ )
1251 {
1252 // Rotation and translation of a plate inside the assembly
1253
1254 G4ThreeVector transl = tgrAssembly->GetComponentPos(ii);
1255 G4String rmName = tgrAssembly->GetComponentRM(ii);
1257 ->FindOrBuildG4RotMatrix( rmName );
1258
1259 //----- Get G4LogicalVolume of component
1260 G4String lvname = tgrAssembly->GetComponentName(ii);
1261 G4LogicalVolume* logvol = g4vmgr->FindG4LogVol( lvname);
1262 if( logvol == 0 )
1263 {
1264 g4vmgr->FindVolume( lvname )->ConstructG4Volumes( 0, 0);
1265 logvol = g4vmgr->FindG4LogVol( lvname, true );
1266 }
1267 // Fill the assembly by the plates
1268 theG4AssemblyVolume->AddPlacedVolume( logvol, transl, rotmat );
1269#ifdef G4VERBOSE
1271 {
1272 G4cout << " G4AssemblyVolume->AddPlacedVolume " << ii
1273 << " " << logvol->GetName()
1274 << " translation " << transl
1275 << " rotmat " << rotmat->colX()
1276 << " " << rotmat->colY()
1277 << " " << rotmat->colZ() << G4endl;
1278 }
1279#endif
1280 }
1281 }
1282
1283 // Rotation and Translation of the assembly inside the world
1284
1285 G4tgrPlaceSimple* placeSimple = (G4tgrPlaceSimple*)place;
1286 G4String rmName = placeSimple->GetRotMatName();
1288 ->FindOrBuildG4RotMatrix( rmName );
1289 G4ThreeVector transl = place->GetPlacement();
1290
1291 G4LogicalVolume* parentLV_nonconst =
1292 const_cast<G4LogicalVolume*>(parentLV);
1293 theG4AssemblyVolume->MakeImprint( parentLV_nonconst, transl, rotmat );
1294
1295 }
1296 else // If it is G4tgrVolumeAssembly
1297 {
1298 G4String ErrMessage = "Volume type not supported: "
1299 + theTgrVolume->GetType() + ", sorry...";
1300 G4Exception("G4tgbVolume::ConstructG4PhysVol()", "NotImplemented",
1301 FatalException, ErrMessage);
1302 }
1303 }
1304
1305 return physvol;
1306}
CLHEP::Hep3Vector G4ThreeVector
HepGeom::Rotate3D G4Rotate3D
G4DLLIMPORT std::ostream G4cerr
G4DivType
@ DivByNdiv
@ DivByWidth
@ DivByNdivAndWidth
Hep3Vector cross(const Hep3Vector &) const
Hep3Vector colX() const
Hep3Vector colY() const
Hep3Vector colZ() const
void MakeImprint(G4LogicalVolume *pMotherLV, G4ThreeVector &translationInMother, G4RotationMatrix *pRotationInMother, G4int copyNumBase=0, G4bool surfCheck=false)
void AddPlacedVolume(G4LogicalVolume *pPlacedVolume, G4ThreeVector &translation, G4RotationMatrix *rotation)
G4VSolid * GetSolid() const
static G4ReflectionFactory * Instance()
G4PhysicalVolumesPair Place(const G4Transform3D &transform3D, const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, G4bool isMany, G4int copyNo, G4bool surfCheck=false)
const G4String & GetName() const
static G4tgbRotationMatrixMgr * GetInstance()
G4RotationMatrix * FindOrBuildG4RotMatrix(const G4String &name)
G4tgbVolume * FindVolume(const G4String &volname)
G4LogicalVolume * FindG4LogVol(const G4String &theName, const G4bool bExists=0)
static G4tgbVolumeMgr * GetInstance()
G4VSolid * BuildSolidForDivision(G4VSolid *parentSolid, EAxis axis)
void ConstructG4Volumes(const G4tgrPlace *place, const G4LogicalVolume *parentLV)
Definition: G4tgbVolume.cc:137
G4double GetOffset() const
G4double GetWidth() const
EAxis GetAxis() const
G4DivType GetDivType() const
G4int GetNDiv() const
const G4String & GetParamType() const
const G4String & GetRotMatName() const
virtual G4ThreeVector GetPlacement() const
Definition: G4tgrPlace.cc:52
const G4String & GetType() const
Definition: G4tgrPlace.hh:61
unsigned int GetCopyNo() const
Definition: G4tgrPlace.hh:60
G4ThreeVector GetComponentPos(G4int ii) const
const G4String & GetComponentRM(G4int ii) const
G4int GetNoComponents() const
const G4String & GetComponentName(G4int ii) const
G4tgrPlaceDivRep * GetPlaceDivision()
G4bool GetCheckOverlaps() const
Definition: G4tgrVolume.hh:100
const G4String & GetType() const
Definition: G4tgrVolume.hh:91
EAxis
Definition: geomdefs.hh:54

Referenced by ConstructG4Volumes().

◆ ConstructG4Volumes()

void G4tgbVolume::ConstructG4Volumes ( const G4tgrPlace place,
const G4LogicalVolume parentLV 
)

Definition at line 137 of file G4tgbVolume.cc.

139{
140#ifdef G4VERBOSE
142 {
143 G4cout << G4endl << "@@@ G4tgbVolume::ConstructG4Volumes - " << GetName() << G4endl;
144 if( place && parentLV ) G4cout << " place in LV " << parentLV->GetName() << G4endl;
145 }
146#endif
148 G4LogicalVolume* logvol = g4vmgr->FindG4LogVol( GetName() );
149 G4bool bFirstCopy = false;
150 if( logvol == 0 )
151 {
152 bFirstCopy = true;
153 if( theTgrVolume->GetType() != "VOLDivision" )
154 {
155 //--- If first time build solid and LogVol
156 G4VSolid* solid = FindOrConstructG4Solid( theTgrVolume->GetSolid() );
157 if( solid != 0 ) // for G4AssemblyVolume it is 0
158 {
159 g4vmgr->RegisterMe( solid );
160 logvol = ConstructG4LogVol( solid );
161 g4vmgr->RegisterMe( logvol );
162 g4vmgr->RegisterChildParentLVs( logvol, parentLV );
163 }
164 else
165 {
166 return;
167 }
168 }
169 else
170 {
171 return;
172 }
173 }
174 //--- Construct PhysVol
175 G4VPhysicalVolume* physvol = ConstructG4PhysVol( place, logvol, parentLV );
176 if( physvol != 0 ) // 0 for G4AssemblyVolumes
177 {
178 g4vmgr->RegisterMe( physvol );
179
180 if( logvol == 0 ) // case of divisions
181 {
182 logvol = physvol->GetLogicalVolume();
183 }
184 }
185 else
186 {
187 return;
188 }
189
190 //--- If first copy build children placements in this LogVol
191 if(bFirstCopy)
192 {
193 std::pair<G4mmapspl::iterator, G4mmapspl::iterator> children
195 G4mmapspl::iterator cite;
196 for( cite = children.first; cite != children.second; cite++ )
197 {
198 //----- Call G4tgrPlace ->constructG4Volumes
199 //---- find G4tgbVolume corresponding to the G4tgrVolume
200 // pointed by G4tgrPlace
201 G4tgrPlace* pl = const_cast<G4tgrPlace*>((*cite).second);
202 G4tgbVolume* svol = g4vmgr->FindVolume( pl->GetVolume()->GetName() );
203 //--- find copyNo
204#ifdef G4VERBOSE
206 {
207 G4cout << " G4tgbVolume::ConstructG4Volumes - construct daughter " << pl->GetVolume()->GetName() << " # " << pl->GetCopyNo() << G4endl;
208 }
209#endif
210 svol->ConstructG4Volumes( pl, logvol );
211 }
212 }
213
214}
bool G4bool
Definition: G4Types.hh:67
G4LogicalVolume * GetLogicalVolume() const
void RegisterMe(const G4tgbVolume *vol)
void RegisterChildParentLVs(const G4LogicalVolume *logvol, const G4LogicalVolume *parentLV)
G4VPhysicalVolume * ConstructG4PhysVol(const G4tgrPlace *place, const G4LogicalVolume *currentLV, const G4LogicalVolume *parentLV)
Definition: G4tgbVolume.cc:942
G4VSolid * FindOrConstructG4Solid(const G4tgrSolid *vol)
Definition: G4tgbVolume.cc:219
G4LogicalVolume * ConstructG4LogVol(const G4VSolid *solid)
Definition: G4tgbVolume.cc:851
G4tgrVolume * GetVolume() const
Definition: G4tgrPlace.hh:59
static G4tgrVolumeMgr * GetInstance()
std::pair< G4mmapspl::iterator, G4mmapspl::iterator > GetChildren(const G4String &name)
const G4String & GetName() const
Definition: G4tgrVolume.hh:89
G4tgrSolid * GetSolid() const
Definition: G4tgrVolume.hh:92

Referenced by G4tgbDetectorConstruction::Construct(), G4tgbDetectorBuilder::ConstructDetector(), ConstructG4PhysVol(), and ConstructG4Volumes().

◆ FindOrConstructG4Solid()

G4VSolid * G4tgbVolume::FindOrConstructG4Solid ( const G4tgrSolid vol)

Definition at line 219 of file G4tgbVolume.cc.

220{
221 G4double angularTolerance = G4GeometryTolerance::GetInstance()
223
224 if( sol == 0 ) { return 0; }
225
226#ifdef G4VERBOSE
228 {
229 G4cout << " G4tgbVolume::FindOrConstructG4Solid():" << G4endl
230 << " SOLID = " << sol << G4endl
231 << " " << sol->GetName() << " of type " << sol->GetType()
232 << G4endl;
233 }
234#endif
235
236 //----- Check if solid exists already
238 ->FindG4Solid( sol->GetName() );
239 if( solid ) { return solid; }
240
241 // Give 'sol' as Boolean solids needs to call this method twice
242
243#ifdef G4VERBOSE
245 {
246 G4cout << " G4tgbVolume::FindOrConstructG4Solid() - "
247 << sol->GetSolidParams().size() << G4endl;
248 }
249#endif
250
251 std::vector<G4double> solParam;
252
253 // In case of BOOLEAN solids, solidParams are taken from components
254
255 if( sol->GetSolidParams().size() == 1)
256 {
257 solParam = * sol->GetSolidParams()[ 0 ];
258 }
259
260 //----------- instantiate the appropiate G4VSolid type
261 G4String stype = sol->GetType();
262 G4String sname = sol->GetName();
263
264 if( stype == "BOX" )
265 {
266 CheckNoSolidParams( stype, 3, solParam.size() );
267 solid = new G4Box( sname, solParam[0], solParam[1], solParam[2] );
268
269 }
270 else if( stype == "TUBE" )
271 {
272 CheckNoSolidParams( stype, 3, solParam.size() );
273 solid = new G4Tubs( sname, solParam[0], solParam[1], solParam[2],
274 0.*deg, 360.*deg );
275 }
276 else if( stype == "TUBS" )
277 {
278 CheckNoSolidParams( stype, 5, solParam.size() );
279 G4double phiDelta = solParam[4];
280 if( std::fabs(phiDelta - twopi) < angularTolerance ) { phiDelta = twopi; }
281 solid = new G4Tubs( sname, solParam[0], solParam[1], solParam[2],
282 solParam[3], phiDelta );
283 }
284 else if( stype == "TRAP" )
285 {
286 if( solParam.size() == 11 )
287 {
288 solid = new G4Trap( sname, solParam[0], solParam[1], solParam[2],
289 solParam[3], solParam[4], solParam[5], solParam[6],
290 solParam[7], solParam[8], solParam[9], solParam[10] );
291 }
292 else if( solParam.size() == 4 )
293 {
294 solid = new G4Trap( sname, solParam[0], solParam[1]/deg,
295 solParam[2]/deg, solParam[3]);
296 }
297 else
298 {
299 G4String ErrMessage1 = "Solid type " + stype;
300 G4String ErrMessage2 = " should have 11 or 4 parameters,\n";
301 G4String ErrMessage3 = "and it has "
302 + G4UIcommand::ConvertToString(G4int(solParam.size()));
303 G4String ErrMessage = ErrMessage1 + ErrMessage2 + ErrMessage3 + " !";
304 G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
305 "InvalidSetup", FatalException, ErrMessage);
306 return 0;
307 }
308
309 }
310 else if( stype == "TRD" )
311 {
312 CheckNoSolidParams( stype, 5, solParam.size() );
313 solid = new G4Trd( sname, solParam[0], solParam[1], solParam[2],
314 solParam[3], solParam[4] );
315 }
316 else if( stype == "PARA" )
317 {
318 CheckNoSolidParams( stype, 6, solParam.size() );
319 solid = new G4Para( sname, solParam[0], solParam[1], solParam[2],
320 solParam[3], solParam[4], solParam[5] );
321 }
322 else if( stype == "CONE" )
323 {
324 CheckNoSolidParams( stype, 5, solParam.size() );
325 solid = new G4Cons( sname, solParam[0], solParam[1], solParam[2],
326 solParam[3], solParam[4], 0., 360.*deg);
327 }
328 else if( stype == "CONS" )
329 {
330 CheckNoSolidParams( stype, 7, solParam.size() );
331 G4double phiDelta = solParam[6];
332 if( std::fabs(phiDelta - twopi) < angularTolerance ) { phiDelta = twopi; }
333 solid = new G4Cons( sname, solParam[0], solParam[1], solParam[2],
334 solParam[3], solParam[4], solParam[5], phiDelta);
335 }
336 else if( stype == "SPHERE" )
337 {
338 CheckNoSolidParams( stype, 6, solParam.size() );
339 G4double phiDelta = solParam[3];
340 if( std::fabs(phiDelta - twopi) < angularTolerance ) { phiDelta = twopi; }
341 G4double thetaDelta = solParam[5];
342 if( std::fabs(thetaDelta - pi) < angularTolerance ) { thetaDelta = pi; }
343 solid = new G4Sphere( sname, solParam[0], solParam[1], solParam[2],
344 phiDelta, solParam[4], thetaDelta);
345 }
346 else if( stype == "ORB" )
347 {
348 CheckNoSolidParams( stype, 1, solParam.size() );
349 solid = new G4Orb( sname, solParam[0] );
350 }
351 else if( stype == "TORUS" )
352 {
353 CheckNoSolidParams( stype, 5, solParam.size() );
354 G4double phiDelta = solParam[4];
355 if( std::fabs(phiDelta - twopi) < angularTolerance ) { phiDelta = twopi; }
356 solid = new G4Torus( sname, solParam[0], solParam[1], solParam[2],
357 solParam[3], phiDelta );
358 }
359 else if( stype == "POLYCONE" )
360 {
361 size_t nplanes = size_t(solParam[2]);
362 G4bool genericPoly = false;
363 if( solParam.size() == 3+nplanes*3 )
364 {
365 genericPoly = true;
366 }
367 else if( solParam.size() == 3+nplanes*2 )
368 {
369 genericPoly = false;
370 }
371 else
372 {
373 G4String Err1 = "Solid type " + stype + " should have ";
374 G4String Err2 = G4UIcommand::ConvertToString(G4int(3+nplanes*3))
375 + " (Z,Rmin,Rmax)\n";
376 G4String Err3 = "or " + G4UIcommand::ConvertToString(G4int(3+nplanes*2));
377 G4String Err4 = " (RZ corners) parameters,\n";
378 G4String Err5 = "and it has "
379 + G4UIcommand::ConvertToString(G4int(solParam.size()));
380 G4String ErrMessage = Err1 + Err2 + Err3 + Err4 + Err5 + " !";
381 G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
382 "InvalidSetup", FatalException, ErrMessage);
383 return 0;
384 }
385
386 if( genericPoly )
387 {
388 std::vector<G4double>* z_p = new std::vector<G4double>;
389 std::vector<G4double>* rmin_p = new std::vector<G4double>;
390 std::vector<G4double>* rmax_p = new std::vector<G4double>;
391 for( size_t ii = 0; ii < nplanes; ii++ )
392 {
393 (*z_p).push_back( solParam[3+3*ii] );
394 (*rmin_p).push_back( solParam[3+3*ii+1] );
395 (*rmax_p).push_back( solParam[3+3*ii+2] );
396 }
397 G4double phiTotal = solParam[1];
398 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
399 solid = new G4Polycone( sname, solParam[0], phiTotal, // start,delta-phi
400 nplanes, // sections
401 &((*z_p)[0]), &((*rmin_p)[0]), &((*rmax_p)[0]));
402 }
403 else
404 {
405 std::vector<G4double>* R_c = new std::vector<G4double>;
406 std::vector<G4double>* Z_c = new std::vector<G4double>;
407 for( size_t ii = 0; ii < nplanes; ii++ )
408 {
409 (*R_c).push_back( solParam[3+2*ii] );
410 (*Z_c).push_back( solParam[3+2*ii+1] );
411 }
412 G4double phiTotal = solParam[1];
413 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
414 solid = new G4Polycone( sname, solParam[0], phiTotal, // start,delta-phi
415 nplanes, // sections
416 &((*R_c)[0]), &((*Z_c)[0]));
417 }
418
419 }
420 else if( stype == "POLYHEDRA" )
421 {
422 size_t nplanes = size_t(solParam[3]);
423 G4bool genericPoly = false;
424 if( solParam.size() == 4+nplanes*3 )
425 {
426 genericPoly = true;
427 }
428 else if( solParam.size() == 4+nplanes*2 )
429 {
430 genericPoly = false;
431 }
432 else
433 {
434 G4String Err1 = "Solid type " + stype + " should have ";
435 G4String Err2 = G4UIcommand::ConvertToString(G4int(4+nplanes*3))
436 + " (Z,Rmin,Rmax)\n";
437 G4String Err3 = "or " + G4UIcommand::ConvertToString(G4int(4+nplanes*2));
438 G4String Err4 = " (RZ corners) parameters,\n";
439 G4String Err5 = "and it has "
440 + G4UIcommand::ConvertToString(G4int(solParam.size()));
441 G4String ErrMessage = Err1 + Err2 + Err3 + Err4 + Err5 + " !";
442 G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
443 "InvalidSetup", FatalException, ErrMessage);
444 return 0;
445 }
446
447 if( genericPoly )
448 {
449 std::vector<G4double>* z_p = new std::vector<G4double>;
450 std::vector<G4double>* rmin_p = new std::vector<G4double>;
451 std::vector<G4double>* rmax_p = new std::vector<G4double>;
452 for( size_t ii = 0; ii < nplanes; ii++ )
453 {
454 (*z_p).push_back( solParam[4+3*ii] );
455 (*rmin_p).push_back( solParam[4+3*ii+1] );
456 (*rmax_p).push_back( solParam[4+3*ii+2] );
457 }
458 G4double phiTotal = solParam[1];
459 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
460 solid = new G4Polyhedra( sname, solParam[0], phiTotal,
461 G4int(solParam[2]), nplanes,
462 &((*z_p)[0]), &((*rmin_p)[0]), &((*rmax_p)[0]));
463 }
464 else
465 {
466 std::vector<G4double>* R_c = new std::vector<G4double>;
467 std::vector<G4double>* Z_c = new std::vector<G4double>;
468 for( size_t ii = 0; ii < nplanes; ii++ )
469 {
470 (*R_c).push_back( solParam[4+2*ii] );
471 (*Z_c).push_back( solParam[4+2*ii+1] );
472 }
473 G4double phiTotal = solParam[1];
474 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
475 solid = new G4Polyhedra( sname, solParam[0], phiTotal,
476 G4int(solParam[2]), nplanes,
477 &((*R_c)[0]), &((*Z_c)[0]));
478 }
479 }
480 else if( stype == "ELLIPTICALTUBE" )
481 {
482 CheckNoSolidParams( stype, 3, solParam.size() );
483 solid = new G4EllipticalTube( sname, solParam[0], solParam[1], solParam[2]);
484 }
485 else if( stype == "ELLIPSOID" )
486 {
487 CheckNoSolidParams( stype, 5, solParam.size() );
488 solid = new G4Ellipsoid( sname, solParam[0], solParam[1], solParam[2],
489 solParam[3], solParam[4] );
490 }
491 else if( stype == "ELLIPTICALCONE" )
492 {
493 CheckNoSolidParams( stype, 4, solParam.size() );
494 solid = new G4EllipticalCone( sname, solParam[0], solParam[1],
495 solParam[2], solParam[3] );
496 }
497 else if( stype == "HYPE" )
498 {
499 CheckNoSolidParams( stype, 5, solParam.size() );
500 solid = new G4Hype( sname, solParam[0], solParam[1], solParam[2],
501 solParam[3], solParam[4] );
502 }
503 else if( stype == "TET" )
504 {
505 CheckNoSolidParams( stype, 12, solParam.size() );
506 G4ThreeVector anchor(solParam[0], solParam[1], solParam[2]);
507 G4ThreeVector p2(solParam[3], solParam[4], solParam[5]);
508 G4ThreeVector p3(solParam[6], solParam[7], solParam[8]);
509 G4ThreeVector p4(solParam[9], solParam[10], solParam[11]);
510 solid = new G4Tet( sname, anchor, p2, p3, p4 );
511 }
512 else if( stype == "TWISTEDBOX" )
513 {
514 CheckNoSolidParams( stype, 4, solParam.size() );
515 solid = new G4TwistedBox( sname, solParam[0], solParam[1],
516 solParam[2], solParam[3]);
517 }
518 else if( stype == "TWISTEDTRAP" )
519 {
520 CheckNoSolidParams( stype, 11, solParam.size() );
521 solid = new G4TwistedTrap( sname, solParam[0], solParam[1], solParam[2],
522 solParam[3], solParam[4], solParam[5], solParam[6],
523 solParam[7], solParam[8], solParam[9], solParam[10] );
524 }
525 else if( stype == "TWISTEDTRD" )
526 {
527 CheckNoSolidParams( stype, 6, solParam.size() );
528 solid = new G4TwistedTrd( sname, solParam[0], solParam[1], solParam[2],
529 solParam[3], solParam[4], solParam[5]);
530 }
531 else if( stype == "TWISTEDTUBS" )
532 {
533 CheckNoSolidParams( stype, 5, solParam.size() );
534 G4double phiTotal = solParam[4];
535 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
536 solid = new G4TwistedTubs( sname, solParam[0], solParam[1], solParam[2],
537 solParam[3], phiTotal);
538 }
539 else if( stype == "BREPBOX" ) // EntityType is = "Closed_Shell"
540 {
541 CheckNoSolidParams( stype, 24, solParam.size() );
542 std::vector<G4Point3D> points;
543 for( size_t ii = 0; ii < 8; ii++ )
544 {
545 points.push_back( G4Point3D(solParam[ii*3+0],
546 solParam[ii*3+1],
547 solParam[ii*3+2]) );
548 }
549 solid = new G4BREPSolidBox( sname, points[0], points[1], points[2],
550 points[3], points[4], points[5], points[6],
551 points[7] );
552 }
553 else if( stype == "BREPCYLINDER" ) // EntityType is = "Closed_Shell"
554 {
555 CheckNoSolidParams( stype, 11, solParam.size() );
556 solid = new G4BREPSolidCylinder( sname,
557 G4ThreeVector( solParam[0], solParam[1], solParam[2] ),
558 G4ThreeVector( solParam[3], solParam[4], solParam[5] ),
559 G4ThreeVector( solParam[6], solParam[7], solParam[8] ),
560 solParam[9], solParam[10] );
561 }
562 else if( stype == "BREPCONE" ) // EntityType is = "Closed_Shell"
563 {
564 CheckNoSolidParams( stype, 12, solParam.size() );
565 solid = new G4BREPSolidCone( sname,
566 G4ThreeVector( solParam[0], solParam[1], solParam[2] ),
567 G4ThreeVector( solParam[3], solParam[4], solParam[5] ),
568 G4ThreeVector( solParam[6], solParam[7], solParam[8] ),
569 solParam[9], solParam[10], solParam[11] );
570 }
571 else if( stype == "BREPSPHERE" ) // EntityType is = "Closed_Shell"
572 {
573 CheckNoSolidParams( stype, 10, solParam.size() );
574 solid = new G4BREPSolidSphere( sname,
575 G4ThreeVector( solParam[0], solParam[1], solParam[2] ),
576 G4ThreeVector( solParam[3], solParam[4], solParam[5] ),
577 G4ThreeVector( solParam[6], solParam[7], solParam[8] ),
578 solParam[9] );
579
580 }
581 else if( stype == "BREPTORUS" ) // EntityType is = "Closed_Shell"
582 {
583 CheckNoSolidParams( stype, 11, solParam.size() );
584 solid = new G4BREPSolidTorus( sname,
585 G4ThreeVector( solParam[0], solParam[1], solParam[2] ),
586 G4ThreeVector( solParam[3], solParam[4], solParam[5] ),
587 G4ThreeVector( solParam[6], solParam[7], solParam[8] ),
588 solParam[9], solParam[10] );
589 }
590 else if( stype == "BREPPCONE" ) // EntityType is = "Closed_Shell"
591 {
592 size_t nplanes = size_t(solParam[2]);
593 CheckNoSolidParams( stype, 4+3*nplanes, solParam.size() );
594 std::vector<G4double>* z_p = new std::vector<G4double>;
595 std::vector<G4double>* rmin_p = new std::vector<G4double>;
596 std::vector<G4double>* rmax_p = new std::vector<G4double>;
597 for( size_t ii = 0; ii < nplanes; ii++ )
598 {
599 (*z_p).push_back( solParam[4+3*ii] );
600 (*rmin_p).push_back( solParam[4+3*ii+1] );
601 (*rmax_p).push_back( solParam[4+3*ii+2] );
602 }
603 G4double phiTotal = solParam[1];
604 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
605 CheckNoSolidParams( stype, 12, solParam.size() );
606 solid = new G4BREPSolidPCone( sname, solParam[0], phiTotal, // start,dph
607 nplanes, // sections
608 solParam[3], // z_start
609 &((*z_p)[0]), &((*rmin_p)[0]),
610 &((*rmax_p)[0]));
611 }
612 else if( stype == "BREPPOLYHEDRA" ) // EntityType is = "Closed_Shell"
613 {
614 size_t nplanes = size_t(solParam[3]);
615 CheckNoSolidParams( stype, 5+3*nplanes, solParam.size() );
616 std::vector<G4double>* z_p = new std::vector<G4double>;
617 std::vector<G4double>* rmin_p = new std::vector<G4double>;
618 std::vector<G4double>* rmax_p = new std::vector<G4double>;
619 for( size_t ii = 0; ii < nplanes; ii++ )
620 {
621 (*z_p).push_back( solParam[5+3*ii] );
622 (*rmin_p).push_back( solParam[5+3*ii+1] );
623 (*rmax_p).push_back( solParam[5+3*ii+2] );
624 }
625 G4double phiTotal = solParam[1];
626 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
627 CheckNoSolidParams( stype, 12, solParam.size() );
628 solid = new G4BREPSolidPolyhedra( sname, solParam[0], phiTotal, // start,dph
629 G4int(solParam[2]), // sides
630 nplanes, // sections
631 solParam[4], // z_start
632 &((*z_p)[0]), &((*rmin_p)[0]),
633 &((*rmax_p)[0]));
634 }
635 else if( stype == "BREPOPENPCONE" ) // EntityType is = "Closed_Shell"
636 {
637 size_t nplanes = size_t(solParam[2]);
638 std::vector<G4double>* z_p = new std::vector<G4double>;
639 std::vector<G4double>* rmin_p = new std::vector<G4double>;
640 std::vector<G4double>* rmax_p = new std::vector<G4double>;
641 for( size_t ii = 0; ii < nplanes; ii++ )
642 {
643 (*z_p).push_back( solParam[4+3*ii] );
644 (*rmin_p).push_back( solParam[4+3*ii+1] );
645 (*rmax_p).push_back( solParam[4+3*ii+2] );
646 }
647 G4double phiTotal = solParam[1];
648 if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
649 CheckNoSolidParams( stype, 12, solParam.size() );
650 solid = new G4BREPSolidOpenPCone( sname, solParam[0], phiTotal, // start,dph
651 nplanes, // sections
652 solParam[3], // z_start
653 &((*z_p)[0]), &((*rmin_p)[0]),
654 &((*rmax_p)[0]));
655 }
656 else if( stype == "TESSELLATED" )
657 {
658 G4int nFacets = G4int(solParam[0]);
659 G4int jj = 0;
660 solid = new G4TessellatedSolid(sname);
661 G4TessellatedSolid* solidTS = (G4TessellatedSolid*)(solid);
662 G4VFacet* facet=0;
663
664 for( G4int ii = 0; ii < nFacets; ii++){
665 G4int nPoints = G4int(solParam[jj+1]);
666 if( G4int(solParam.size()) < jj + nPoints*3 + 2 ) {
667 G4String Err1 = "Too small number of parameters in tesselated solid, it should be at least " + G4UIcommand::ConvertToString(jj + nPoints*3 + 2 );
668 G4String Err2 = " facet number " + G4UIcommand::ConvertToString(ii);
669 G4String Err3 = " number of parameters is " + G4UIcommand::ConvertToString(G4int(solParam.size()));
670 G4String ErrMessage = Err1 + Err2 + Err3 + " !";
671 G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
672 "InvalidSetup", FatalException, ErrMessage);
673 return 0;
674 }
675
676 if( nPoints == 3 )
677 {
678 G4ThreeVector pt0(solParam[jj+2],solParam[jj+3],solParam[jj+4]);
679 G4ThreeVector vt1(solParam[jj+5],solParam[jj+6],solParam[jj+7]);
680 G4ThreeVector vt2(solParam[jj+8],solParam[jj+9],solParam[jj+10]);
681 G4FacetVertexType vertexType = ABSOLUTE;
682 if( solParam[jj+11] == 0 )
683 {
684 vertexType = ABSOLUTE;
685 }
686 else if( solParam[jj+11] == 1 )
687 {
688 vertexType = RELATIVE;
689 }
690 else
691 {
692 G4String Err1 = "Wrong number of vertex type in tesselated solid, it should be 0 =ABSOLUTE) or 1 (=RELATIVE)";
693 G4String Err2 = " facet number " + G4UIcommand::ConvertToString(G4int(ii));
694 G4String Err3 = " vertex type is " + G4UIcommand::ConvertToString(solParam[jj+11]);
695 G4String ErrMessage = Err1 + Err2 + Err3 + " !";
696 G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
697 "InvalidSetup", FatalException, ErrMessage);
698 return 0;
699 }
700 facet = new G4TriangularFacet( pt0, vt1, vt2, vertexType );
701 }
702 else if( nPoints == 4 )
703 {
704 G4ThreeVector pt0(solParam[jj+2],solParam[jj+3],solParam[jj+4]);
705 G4ThreeVector vt1(solParam[jj+5],solParam[jj+6],solParam[jj+7]);
706 G4ThreeVector vt2(solParam[jj+8],solParam[jj+9],solParam[jj+10]);
707 G4ThreeVector vt3(solParam[jj+11],solParam[jj+12],solParam[jj+13]);
708 G4FacetVertexType vertexType = ABSOLUTE;
709 if( solParam[jj+14] == 0 )
710 {
711 vertexType = ABSOLUTE;
712 }
713 else if( solParam[jj+14] == 1 )
714 {
715 vertexType = RELATIVE;
716 }
717 else
718 {
719 G4String Err1 = "Wrong number of vertex type in tesselated solid, it should be 0 =ABSOLUTE) or 1 (=RELATIVE)";
720 G4String Err2 = " facet number " + G4UIcommand::ConvertToString(G4int(ii));
721 G4String Err3 = " vertex type is " + G4UIcommand::ConvertToString(solParam[jj+14]);
722 G4String ErrMessage = Err1 + Err2 + Err3 + " !";
723 G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
724 "InvalidSetup", FatalException, ErrMessage);
725 return 0;
726 }
727 facet = new G4QuadrangularFacet( pt0, vt1, vt2, vt3, vertexType );
728 }
729 else
730 {
731 G4String Err1 = "Wrong number of points in tesselated solid, it should be 3 or 4";
732 G4String Err2 = " facet number " + G4UIcommand::ConvertToString(G4int(ii));
733 G4String Err3 = " number of points is " + G4UIcommand::ConvertToString(G4int(nPoints));
734 G4String ErrMessage = Err1 + Err2 + Err3 + " !";
735 G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
736 "InvalidSetup", FatalException, ErrMessage);
737 return 0;
738 }
739
740 solidTS->AddFacet( facet );
741 jj += nPoints*3 + 2;
742 }
743
744 }
745 else if( stype == "EXTRUDED" )
746 {
747 std::vector<G4TwoVector> polygonList;
748 std::vector<G4ExtrudedSolid::ZSection> zsectionList;
749 G4int nPolygons = G4int(solParam[0]);
750 G4int ii = 1;
751 G4int nMax = nPolygons*2+1;
752 for( ;ii < nMax; ii+=2 )
753 {
754 polygonList.push_back( G4TwoVector(solParam[ii],solParam[ii+1]) );
755 }
756 G4int nZSections = G4int(solParam[ii]);
757 nMax = nPolygons*2 + nZSections*4 + 2;
758 ii++;
759 for( ; ii < nMax; ii+=4 )
760 {
761 G4TwoVector offset(solParam[ii+1],solParam[ii+2]);
762 zsectionList.push_back( G4ExtrudedSolid::ZSection(solParam[ii],offset,solParam[ii+3]) );
763 }
764 solid = new G4ExtrudedSolid( sname, polygonList, zsectionList );
765
766 }
767 else if( stype.substr(0,7) == "Boolean" )
768 {
769 const G4tgrSolidBoolean* solb = dynamic_cast<const G4tgrSolidBoolean*>(sol);
770 if (!solb)
771 {
772 G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
773 "InvalidSetup", FatalException, "Invalid Solid pointer");
774 return 0;
775 }
776 G4VSolid* sol1 = FindOrConstructG4Solid( solb->GetSolid(0));
777 G4VSolid* sol2 = FindOrConstructG4Solid( solb->GetSolid(1));
779 ->FindOrBuildG4RotMatrix( sol->GetRelativeRotMatName() );
780 G4ThreeVector relPlace = solb->GetRelativePlace();
781
782 if( stype == "Boolean_UNION" )
783 {
784 solid = new G4UnionSolid( sname, sol1, sol2, relRotMat, relPlace );
785 }
786 else if( stype == "Boolean_SUBTRACTION" )
787 {
788 solid = new G4SubtractionSolid( sname, sol1, sol2, relRotMat, relPlace );
789 }
790 else if( stype == "Boolean_INTERSECTION" )
791 {
792 solid = new G4IntersectionSolid( sname, sol1, sol2, relRotMat, relPlace );
793 }
794 else
795 {
796 G4String ErrMessage = "Unknown Boolean type " + stype;
797 G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
798 "InvalidSetup", FatalException, ErrMessage);
799 return 0;
800 }
801 }
802 else
803 {
804 G4String ErrMessage = "Solids of type " + stype
805 + " not implemented yet, sorry...";
806 G4Exception("G4tgbVolume::FindOrConstructG4Solid()", "NotImplemented",
807 FatalException, ErrMessage);
808 return 0;
809 }
810
811#ifdef G4VERBOSE
813 {
814 G4cout << " G4tgbVolume::FindOrConstructG4Solid()" << G4endl
815 << " Created solid " << sname
816 << " of type " << solid->GetEntityType() << G4endl;
817 }
818#endif
819
820#ifdef G4VERBOSE
822 {
823 G4cout << " Constructing new G4Solid: "
824 << *solid << G4endl;
825 }
826#endif
827
828 return solid;
829}
HepGeom::Point3D< G4double > G4Point3D
Definition: G4Point3D.hh:35
CLHEP::Hep2Vector G4TwoVector
Definition: G4TwoVector.hh:42
G4FacetVertexType
Definition: G4VFacet.hh:56
@ ABSOLUTE
Definition: G4VFacet.hh:56
@ RELATIVE
Definition: G4VFacet.hh:56
G4double GetAngularTolerance() const
Definition: G4Hype.hh:67
Definition: G4Orb.hh:52
G4bool AddFacet(G4VFacet *aFacet)
Definition: G4Tet.hh:57
G4VSolid * FindG4Solid(const G4String &name)
void CheckNoSolidParams(const G4String &solidType, const unsigned int NoParamExpected, const unsigned int NoParam)
Definition: G4tgbVolume.cc:832
const G4tgrSolid * GetSolid(G4int ii) const
G4ThreeVector GetRelativePlace() const
const G4double pi

Referenced by ConstructG4Volumes(), and FindOrConstructG4Solid().

◆ GetColour()

const G4double * G4tgbVolume::GetColour ( ) const
inline

Definition at line 110 of file G4tgbVolume.hh.

110{ return theTgrVolume->GetColour(); }
G4double * GetColour() const
Definition: G4tgrVolume.hh:97

Referenced by ConstructG4LogVol().

◆ GetName()

const G4String & G4tgbVolume::GetName ( ) const
inline

Definition at line 108 of file G4tgbVolume.hh.

108{ return theTgrVolume->GetName(); }

Referenced by BuildSolidForDivision(), ConstructG4LogVol(), ConstructG4PhysVol(), ConstructG4Volumes(), and G4tgbVolumeMgr::RegisterMe().

◆ GetVisibility()

G4bool G4tgbVolume::GetVisibility ( ) const
inline

Definition at line 109 of file G4tgbVolume.hh.

109{ return theTgrVolume->GetVisibility(); }
G4bool GetVisibility() const
Definition: G4tgrVolume.hh:96

Referenced by ConstructG4LogVol().

◆ SetCutsInEnergy()

void G4tgbVolume::SetCutsInEnergy ( G4LogicalVolume logvol,
std::map< G4String, G4double cuts 
)

◆ SetCutsInRange()

void G4tgbVolume::SetCutsInRange ( G4LogicalVolume logvol,
std::map< G4String, G4double cuts 
)

The documentation for this class was generated from the following files: