Geant4 10.7.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 58 of file G4tgbVolume.hh.

Constructor & Destructor Documentation

◆ G4tgbVolume() [1/2]

G4tgbVolume::G4tgbVolume ( )

Definition at line 105 of file G4tgbVolume.cc.

106{
107}

◆ ~G4tgbVolume()

G4tgbVolume::~G4tgbVolume ( )

Definition at line 110 of file G4tgbVolume.cc.

111{
112}

◆ G4tgbVolume() [2/2]

G4tgbVolume::G4tgbVolume ( G4tgrVolume vol)

Definition at line 115 of file G4tgbVolume.cc.

116{
117 theTgrVolume = vol;
118}

Member Function Documentation

◆ BuildSolidForDivision()

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

Definition at line 1207 of file G4tgbVolume.cc.

1208{
1209 G4VSolid* solid = nullptr;
1210 G4double redf =
1211 (parentSolid->GetExtent().GetXmax() - parentSolid->GetExtent().GetXmin());
1212 redf = std::min(redf, parentSolid->GetExtent().GetYmax() -
1213 parentSolid->GetExtent().GetYmin());
1214 redf = std::min(redf, parentSolid->GetExtent().GetZmax() -
1215 parentSolid->GetExtent().GetZmin());
1216 redf *= 0.001; // make daugther much smaller, to fit in parent
1217
1218 if(parentSolid->GetEntityType() == "G4Box")
1219 {
1220 G4Box* psolid = (G4Box*) (parentSolid);
1221 solid = new G4Box(GetName(), psolid->GetXHalfLength() * redf,
1222 psolid->GetZHalfLength() * redf,
1223 psolid->GetZHalfLength() * redf);
1224 }
1225 else if(parentSolid->GetEntityType() == "G4Tubs")
1226 {
1227 G4Tubs* psolid = (G4Tubs*) (parentSolid);
1228 solid = new G4Tubs(GetName(), psolid->GetInnerRadius() * redf,
1229 psolid->GetOuterRadius() * redf,
1230 psolid->GetZHalfLength() * redf,
1231 psolid->GetStartPhiAngle(), psolid->GetDeltaPhiAngle());
1232 }
1233 else if(parentSolid->GetEntityType() == "G4Cons")
1234 {
1235 G4Cons* psolid = (G4Cons*) (parentSolid);
1236 solid = new G4Cons(GetName(), psolid->GetInnerRadiusMinusZ() * redf,
1237 psolid->GetOuterRadiusMinusZ() * redf,
1238 psolid->GetInnerRadiusPlusZ() * redf,
1239 psolid->GetOuterRadiusPlusZ() * redf,
1240 psolid->GetZHalfLength() * redf,
1241 psolid->GetStartPhiAngle(), psolid->GetDeltaPhiAngle());
1242 }
1243 else if(parentSolid->GetEntityType() == "G4Trd")
1244 {
1245 G4Trd* psolid = (G4Trd*) (parentSolid);
1246 G4double mpDx1 = psolid->GetXHalfLength1();
1247 G4double mpDx2 = psolid->GetXHalfLength2();
1248
1249 if(axis == kXAxis &&
1250 std::fabs(mpDx1 - mpDx2) >
1252 {
1253 solid = new G4Trap(GetName(), psolid->GetZHalfLength() * redf,
1254 psolid->GetYHalfLength1() * redf,
1255 psolid->GetXHalfLength2() * redf,
1256 psolid->GetXHalfLength1() * redf);
1257 }
1258 else
1259 {
1260 solid = new G4Trd(
1261 GetName(), psolid->GetXHalfLength1() * redf,
1262 psolid->GetXHalfLength2() * redf, psolid->GetYHalfLength1() * redf,
1263 psolid->GetYHalfLength2() * redf, psolid->GetZHalfLength() * redf);
1264 }
1265 }
1266 else if(parentSolid->GetEntityType() == "G4Para")
1267 {
1268 G4Para* psolid = (G4Para*) (parentSolid);
1269 solid = new G4Para(
1270 GetName(), psolid->GetXHalfLength() * redf,
1271 psolid->GetYHalfLength() * redf, psolid->GetZHalfLength() * redf,
1272 std::atan(psolid->GetTanAlpha()), psolid->GetSymAxis().theta(),
1273 psolid->GetSymAxis().phi());
1274 }
1275 else if(parentSolid->GetEntityType() == "G4Polycone")
1276 {
1277 G4Polycone* psolid = (G4Polycone*) (parentSolid);
1278 G4PolyconeHistorical origParam = *(psolid->GetOriginalParameters());
1279 for(G4int ii = 0; ii < origParam.Num_z_planes; ++ii)
1280 {
1281 origParam.Rmin[ii] = origParam.Rmin[ii] * redf;
1282 origParam.Rmax[ii] = origParam.Rmax[ii] * redf;
1283 }
1284 solid = new G4Polycone(GetName(), psolid->GetStartPhi(),
1285 psolid->GetEndPhi(), origParam.Num_z_planes,
1286 origParam.Z_values, origParam.Rmin, origParam.Rmax);
1287 }
1288 else if(parentSolid->GetEntityType() == "G4GenericPolycone")
1289 {
1290 G4GenericPolycone* psolid = (G4GenericPolycone*) (parentSolid);
1291 const G4int numRZ = psolid->GetNumRZCorner();
1292 G4double* r = new G4double[numRZ];
1293 G4double* z = new G4double[numRZ];
1294 for(G4int ii = 0; ii < numRZ; ++ii)
1295 {
1296 r[ii] = psolid->GetCorner(ii).r;
1297 z[ii] = psolid->GetCorner(ii).z;
1298 }
1299 solid = new G4GenericPolycone(GetName(), psolid->GetStartPhi(),
1300 psolid->GetEndPhi() - psolid->GetStartPhi(),
1301 numRZ, r, z);
1302 delete[] r;
1303 delete[] z;
1304 }
1305 else if(parentSolid->GetEntityType() == "G4Polyhedra")
1306 {
1307 G4Polyhedra* psolid = (G4Polyhedra*) (parentSolid);
1308 G4PolyhedraHistorical origParam = *(psolid->GetOriginalParameters());
1309 for(G4int ii = 0; ii < origParam.Num_z_planes; ++ii)
1310 {
1311 origParam.Rmin[ii] = origParam.Rmin[ii] * redf;
1312 origParam.Rmax[ii] = origParam.Rmax[ii] * redf;
1313 }
1314 solid =
1315 new G4Polyhedra(GetName(), psolid->GetStartPhi(), psolid->GetEndPhi(),
1316 psolid->GetNumSide(), origParam.Num_z_planes,
1317 origParam.Z_values, origParam.Rmin, origParam.Rmax);
1318 }
1319 else
1320 {
1321 G4String ErrMessage = "Solid type not supported. VOLUME= " + GetName() +
1322 " Solid type= " + parentSolid->GetEntityType() +
1323 "\n" +
1324 "Only supported types are: G4Box, G4Tubs, G4Cons," +
1325 " G4Trd, G4Para, G4Polycone, G4Polyhedra.";
1326 G4Exception("G4tgbVolume::BuildSolidForDivision()", "NotImplemented",
1327 FatalException, ErrMessage);
1328 return nullptr;
1329 }
1330
1331#ifdef G4VERBOSE
1333 {
1334 G4cout << " Constructing new G4Solid for division: " << *solid << G4endl;
1335 }
1336#endif
1337 return solid;
1338}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double phi() const
double theta() const
Definition: G4Box.hh:56
G4double GetZHalfLength() const
G4double GetXHalfLength() const
Definition: G4Cons.hh:78
G4double GetOuterRadiusPlusZ() const
G4double GetStartPhiAngle() const
G4double GetDeltaPhiAngle() const
G4double GetInnerRadiusMinusZ() const
G4double GetInnerRadiusPlusZ() const
G4double GetOuterRadiusMinusZ() const
G4double GetZHalfLength() const
G4double GetStartPhi() const
G4double GetEndPhi() const
G4int GetNumRZCorner() const
G4PolyconeSideRZ GetCorner(G4int index) const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
Definition: G4Para.hh:79
G4double GetTanAlpha() const
G4ThreeVector GetSymAxis() const
G4double GetYHalfLength() const
G4double GetZHalfLength() const
G4double GetXHalfLength() const
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:75
G4double GetZHalfLength() const
G4double GetInnerRadius() const
G4double GetOuterRadius() const
G4double GetStartPhiAngle() const
G4double GetDeltaPhiAngle() const
virtual G4VisExtent GetExtent() const
Definition: G4VSolid.cc:670
virtual G4GeometryType GetEntityType() const =0
G4double GetYmin() const
Definition: G4VisExtent.hh:101
G4double GetXmax() const
Definition: G4VisExtent.hh:100
G4double GetYmax() const
Definition: G4VisExtent.hh:102
G4double GetZmax() const
Definition: G4VisExtent.hh:104
G4double GetZmin() const
Definition: G4VisExtent.hh:103
G4double GetXmin() const
Definition: G4VisExtent.hh:99
const G4String & GetName() const
Definition: G4tgbVolume.hh:99
static G4int GetVerboseLevel()
@ kXAxis
Definition: geomdefs.hh:55

Referenced by ConstructG4PhysVol().

◆ CheckNoSolidParams()

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

Definition at line 747 of file G4tgbVolume.cc.

750{
751 if(NoParamExpected != NoParam)
752 {
753 G4String Err1 = "Solid type " + solidType + " should have ";
754 G4String Err2 =
755 G4UIcommand::ConvertToString(G4int(NoParamExpected)) + " parameters,\n";
756 G4String Err3 =
757 "and it has " + G4UIcommand::ConvertToString(G4int(NoParam));
758 G4String ErrMessage = Err1 + Err2 + Err3 + " !";
759 G4Exception("G4tgbVolume::CheckNoSolidParams()", "InvalidSetup",
760 FatalException, ErrMessage);
761 }
762}
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:430

Referenced by FindOrConstructG4Solid().

◆ ConstructG4LogVol()

G4LogicalVolume * G4tgbVolume::ConstructG4LogVol ( const G4VSolid solid)

Definition at line 765 of file G4tgbVolume.cc.

766{
767 G4LogicalVolume* logvol;
768
769#ifdef G4VERBOSE
771 {
772 G4cout << " G4tgbVolume::ConstructG4LogVol() - " << GetName() << G4endl;
773 }
774#endif
775
776 //----------- Get the material first
778 theTgrVolume->GetMaterialName());
779 if(mate == nullptr)
780 {
781 G4String ErrMessage = "Material not found " +
782 theTgrVolume->GetMaterialName() + " for volume " +
783 GetName() + ".";
784 G4Exception("G4tgbVolume::ConstructG4LogVol()", "InvalidSetup",
785 FatalException, ErrMessage);
786 }
787#ifdef G4VERBOSE
789 {
790 G4cout << " G4tgbVolume::ConstructG4LogVol() -"
791 << " Material constructed: " << mate->GetName() << G4endl;
792 }
793#endif
794
795 //---------- Construct the LV
796 logvol = new G4LogicalVolume(const_cast<G4VSolid*>(solid),
797 const_cast<G4Material*>(mate), GetName());
798
799#ifdef G4VERBOSE
801 {
802 G4cout << " Constructing new G4LogicalVolume: " << logvol->GetName()
803 << " mate " << mate->GetName() << G4endl;
804 }
805#endif
806
807 //---------- Set visibility and colour
808 if(!GetVisibility() || GetColour()[0] != -1)
809 {
810 G4VisAttributes* visAtt = new G4VisAttributes();
811#ifdef G4VERBOSE
813 {
814 G4cout << " Constructing new G4VisAttributes: " << *visAtt << G4endl;
815 }
816#endif
817
818 if(!GetVisibility())
819 {
820 visAtt->SetVisibility(GetVisibility());
821 }
822 else if(GetColour()[0] != -1)
823 {
824 // this else should not be necessary, because if the visibility
825 // is set to off, colour should have no effect. But it does not
826 // work: if you set colour and vis off, it is visualized!?!?!?
827
828 const G4double* col = GetColour();
829 if(col[3] == -1.)
830 {
831 visAtt->SetColour(G4Colour(col[0], col[1], col[2]));
832 }
833 else
834 {
835 visAtt->SetColour(G4Colour(col[0], col[1], col[2], col[3]));
836 }
837 }
838 logvol->SetVisAttributes(visAtt);
839 }
840
841#ifdef G4VERBOSE
843 {
844 G4cout << " G4tgbVolume::ConstructG4LogVol() -"
845 << " Created logical volume: " << GetName() << G4endl;
846 }
847#endif
848
849 return logvol;
850}
void SetVisAttributes(const G4VisAttributes *pVA)
const G4String & GetName() const
const G4String & GetName() const
Definition: G4Material.hh:175
void SetColour(const G4Colour &)
void SetVisibility(G4bool=true)
G4Material * FindOrBuildG4Material(const G4String &name, G4bool bMustExist=true)
static G4tgbMaterialMgr * GetInstance()
const G4double * GetColour() const
Definition: G4tgbVolume.hh:101
G4bool GetVisibility() const
Definition: G4tgbVolume.hh:100
const G4String & GetMaterialName() const
Definition: G4tgrVolume.hh:87

Referenced by ConstructG4Volumes().

◆ ConstructG4PhysVol()

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

Definition at line 854 of file G4tgbVolume.cc.

857{
858 G4VPhysicalVolume* physvol = nullptr;
859 G4int copyNo;
860
861 //----- Case of placement of top volume
862 if(place == nullptr)
863 {
864#ifdef G4VERBOSE
866 {
867 G4cout << " G4tgbVolume::ConstructG4PhysVol() - World: " << GetName()
868 << G4endl;
869 }
870#endif
871 physvol = new G4PVPlacement(
872 nullptr, G4ThreeVector(), const_cast<G4LogicalVolume*>(currentLV),
873 GetName(), 0, false, 0, theTgrVolume->GetCheckOverlaps());
874#ifdef G4VERBOSE
876 {
877 G4cout << " Constructing new : G4PVPlacement " << physvol->GetName()
878 << G4endl;
879 }
880#endif
881 }
882 else
883 {
884 copyNo = place->GetCopyNo();
885
886#ifdef G4VERBOSE
888 {
889 G4cout << " G4tgbVolume::ConstructG4PhysVol() - " << GetName() << G4endl
890 << " inside " << parentLV->GetName() << " copy No: " << copyNo
891 << " of type= " << theTgrVolume->GetType() << G4endl
892 << " placement type= " << place->GetType() << G4endl;
893 }
894#endif
895
896 if(theTgrVolume->GetType() == "VOLSimple")
897 {
898 //----- Get placement
899#ifdef G4VERBOSE
901 {
902 G4cout << " G4tgbVolume::ConstructG4PhysVol() - Placement type = "
903 << place->GetType() << G4endl;
904 }
905#endif
906
907 //--------------- If it is G4tgrPlaceSimple
908 if(place->GetType() == "PlaceSimple")
909 {
910 //----- Get rotation matrix
911 G4tgrPlaceSimple* placeSimple = (G4tgrPlaceSimple*) place;
912 G4String rmName = placeSimple->GetRotMatName();
913
914 G4RotationMatrix* rotmat =
916 //----- Place volume in mother
917 G4double check =
918 (rotmat->colX().cross(rotmat->colY())) * rotmat->colZ();
919 G4double tol = 1.0e-3;
920 //---- Check that matrix is ortogonal
921 if(1 - std::abs(check) > tol)
922 {
923 G4cerr << " Matrix : " << rmName << " " << rotmat->colX() << " "
924 << rotmat->colY() << " " << rotmat->colZ() << G4endl
925 << " product x X y * z = " << check << " x X y "
926 << rotmat->colX().cross(rotmat->colY()) << G4endl;
927 G4String ErrMessage = "Rotation is not ortogonal " + rmName + " !";
928 G4Exception("G4tgbVolume::ConstructG4PhysVol()", "InvalidSetup",
929 FatalException, ErrMessage);
930 //---- Check if it is reflection
931 }
932 else if(1 + check <= tol)
933 {
934 G4Translate3D transl = place->GetPlacement();
935 G4Transform3D trfrm = transl * G4Rotate3D(*rotmat);
936 physvol =
938 trfrm, GetName(), const_cast<G4LogicalVolume*>(currentLV),
939 const_cast<G4LogicalVolume*>(parentLV), false, copyNo, false))
940 .first;
941 }
942 else
943 {
944#ifdef G4VERBOSE
946 {
947 G4cout << "Construction new G4VPhysicalVolume"
948 << " through G4ReflectionFactory " << GetName()
949 << " in volume " << parentLV->GetName() << " copyNo "
950 << copyNo << " position " << place->GetPlacement() << " ROT "
951 << rotmat->colX() << " " << rotmat->colY() << " "
952 << rotmat->colZ() << G4endl;
953 }
954#endif
955 physvol =
956 new G4PVPlacement(rotmat, place->GetPlacement(),
957 const_cast<G4LogicalVolume*>(currentLV),
958 GetName(), const_cast<G4LogicalVolume*>(parentLV),
959 false, copyNo, theTgrVolume->GetCheckOverlaps());
960 }
961
962 //--------------- If it is G4tgrPlaceParam
963 }
964 else if(place->GetType() == "PlaceParam")
965 {
967
968 //----- See what parameterisation type
969#ifdef G4VERBOSE
971 {
972 G4cout << " G4tgbVolume::ConstructG4PhysVol() -" << G4endl
973 << " param: " << GetName() << " in " << parentLV->GetName()
974 << " param type= " << dp->GetParamType() << G4endl;
975 }
976#endif
977
978 G4tgbPlaceParameterisation* param = nullptr;
979
980 if((dp->GetParamType() == "CIRCLE") ||
981 (dp->GetParamType() == "CIRCLE_XY") ||
982 (dp->GetParamType() == "CIRCLE_XZ") ||
983 (dp->GetParamType() == "CIRCLE_YZ"))
984 {
985 param = new G4tgbPlaceParamCircle(dp);
986 }
987 else if((dp->GetParamType() == "LINEAR") ||
988 (dp->GetParamType() == "LINEAR_X") ||
989 (dp->GetParamType() == "LINEAR_Y") ||
990 (dp->GetParamType() == "LINEAR_Z"))
991 {
992 param = new G4tgbPlaceParamLinear(dp);
993 }
994 else if((dp->GetParamType() == "SQUARE") ||
995 (dp->GetParamType() == "SQUARE_XY") ||
996 (dp->GetParamType() == "SQUARE_XZ") ||
997 (dp->GetParamType() == "SQUARE_YZ"))
998 {
999 param = new G4tgbPlaceParamSquare(dp);
1000 }
1001 else
1002 {
1003 G4String ErrMessage = "Parameterisation has wrong type, TYPE: " +
1004 G4String(dp->GetParamType()) + " !";
1005 G4Exception("G4tgbVolume::ConstructG4PhysVol", "WrongArgument",
1006 FatalException, ErrMessage);
1007 return nullptr;
1008 }
1009#ifdef G4VERBOSE
1011 {
1012 G4cout << " G4tgbVolume::ConstructG4PhysVol() -" << G4endl
1013 << " New G4PVParameterised: " << GetName() << " vol "
1014 << currentLV->GetName() << " in vol " << parentLV->GetName()
1015 << " axis " << param->GetAxis() << " nCopies "
1016 << param->GetNCopies() << G4endl;
1017 }
1018#endif
1019 physvol = new G4PVParameterised(
1020 GetName(), const_cast<G4LogicalVolume*>(currentLV),
1021 const_cast<G4LogicalVolume*>(parentLV), EAxis(param->GetAxis()),
1022 param->GetNCopies(), param);
1023#ifdef G4VERBOSE
1025 {
1026 G4cout << " Constructing new G4PVParameterised: "
1027 << physvol->GetName() << " in volume " << parentLV->GetName()
1028 << " N copies " << param->GetNCopies() << " axis "
1029 << param->GetAxis() << G4endl;
1030 }
1031#endif
1032 }
1033 else if(place->GetType() == "PlaceReplica")
1034 {
1035 //--------------- If it is PlaceReplica
1036 G4tgrPlaceDivRep* dpr = (G4tgrPlaceDivRep*) place;
1037
1038#ifdef G4VERBOSE
1040 {
1041 G4cout << " G4tgbVolume::ConstructG4PhysVol() -" << G4endl
1042 << " replica"
1043 << " " << currentLV->GetName() << " in " << parentLV->GetName()
1044 << " NDiv " << dpr->GetNDiv() << " Width " << dpr->GetWidth()
1045 << " offset " << dpr->GetOffset() << G4endl;
1046 }
1047#endif
1048 physvol = new G4PVReplica(
1049 GetName(), const_cast<G4LogicalVolume*>(currentLV),
1050 const_cast<G4LogicalVolume*>(parentLV), EAxis(dpr->GetAxis()),
1051 dpr->GetNDiv(), dpr->GetWidth(), dpr->GetOffset());
1052#ifdef G4VERBOSE
1054 {
1055 G4cout << " Constructing new G4PVReplica: " << currentLV->GetName()
1056 << " in " << parentLV->GetName() << " NDiv " << dpr->GetNDiv()
1057 << " Width " << dpr->GetWidth() << " offset "
1058 << dpr->GetOffset() << G4endl;
1059 }
1060#endif
1061 }
1062 }
1063 else if(theTgrVolume->GetType() == "VOLDivision")
1064 {
1065 G4tgrVolumeDivision* volr = (G4tgrVolumeDivision*) theTgrVolume;
1066 G4tgrPlaceDivRep* placeDiv = volr->GetPlaceDivision();
1067 G4VSolid* solid =
1068 BuildSolidForDivision(parentLV->GetSolid(), placeDiv->GetAxis());
1070 theTgrVolume->GetMaterialName());
1071 G4LogicalVolume* divLV =
1072 new G4LogicalVolume(solid, const_cast<G4Material*>(mate), GetName());
1073#ifdef G4VERBOSE
1075 {
1076 G4cout << " Constructed new G4LogicalVolume for division: "
1077 << divLV->GetName() << " mate " << mate->GetName() << G4endl;
1078 }
1079#endif
1080
1081 G4DivType divType = placeDiv->GetDivType();
1082 switch(divType)
1083 {
1084 case DivByNdiv:
1085 physvol = new G4PVDivision(GetName(), (G4LogicalVolume*) divLV,
1086 const_cast<G4LogicalVolume*>(parentLV),
1087 placeDiv->GetAxis(), placeDiv->GetNDiv(),
1088 placeDiv->GetOffset());
1089#ifdef G4VERBOSE
1091 {
1092 G4cout << " Constructing new G4PVDivision by number of divisions: "
1093 << GetName() << " in " << parentLV->GetName() << " axis "
1094 << placeDiv->GetAxis() << " Ndiv " << placeDiv->GetNDiv()
1095 << " offset " << placeDiv->GetOffset() << G4endl;
1096 }
1097#endif
1098 break;
1099 case DivByWidth:
1100 physvol = new G4PVDivision(GetName(), (G4LogicalVolume*) divLV,
1101 const_cast<G4LogicalVolume*>(parentLV),
1102 placeDiv->GetAxis(), placeDiv->GetWidth(),
1103 placeDiv->GetOffset());
1104#ifdef G4VERBOSE
1106 {
1107 G4cout << " Constructing new G4PVDivision by width: " << GetName()
1108 << " in " << parentLV->GetName() << " axis "
1109 << placeDiv->GetAxis() << " width " << placeDiv->GetWidth()
1110 << " offset " << placeDiv->GetOffset() << G4endl;
1111 }
1112#endif
1113 break;
1114 case DivByNdivAndWidth:
1115 physvol = new G4PVDivision(
1116 GetName(), (G4LogicalVolume*) divLV,
1117 const_cast<G4LogicalVolume*>(parentLV), placeDiv->GetAxis(),
1118 placeDiv->GetNDiv(), placeDiv->GetWidth(), placeDiv->GetOffset());
1119#ifdef G4VERBOSE
1121 {
1122 G4cout << " Constructing new G4PVDivision by width"
1123 << " and number of divisions: " << GetName() << " in "
1124 << parentLV->GetName() << " axis " << placeDiv->GetAxis()
1125 << " Ndiv " << placeDiv->GetNDiv() << " width "
1126 << placeDiv->GetWidth() << " offset "
1127 << placeDiv->GetOffset() << G4endl;
1128 }
1129#endif
1130 break;
1131 }
1132 }
1133 else if(theTgrVolume->GetType() == "VOLAssembly")
1134 {
1135 // Define one layer as one assembly volume
1136 G4tgrVolumeAssembly* tgrAssembly = (G4tgrVolumeAssembly*) theTgrVolume;
1137
1138 if(!theG4AssemblyVolume)
1139 {
1140 theG4AssemblyVolume = new G4AssemblyVolume;
1141#ifdef G4VERBOSE
1143 {
1144 G4cout << " Constructing new G4AssemblyVolume: "
1145 << " number of assembly components "
1146 << tgrAssembly->GetNoComponents() << G4endl;
1147 }
1148#endif
1150 for(G4int ii = 0; ii < tgrAssembly->GetNoComponents(); ++ii)
1151 {
1152 // Rotation and translation of a plate inside the assembly
1153
1154 G4ThreeVector transl = tgrAssembly->GetComponentPos(ii);
1155 G4String rmName = tgrAssembly->GetComponentRM(ii);
1156 G4RotationMatrix* rotmat =
1158 rmName);
1159
1160 //----- Get G4LogicalVolume of component
1161 G4String lvname = tgrAssembly->GetComponentName(ii);
1162 G4LogicalVolume* logvol = g4vmgr->FindG4LogVol(lvname);
1163 if(logvol == nullptr)
1164 {
1165 g4vmgr->FindVolume(lvname)->ConstructG4Volumes(0, 0);
1166 logvol = g4vmgr->FindG4LogVol(lvname, true);
1167 }
1168 // Fill the assembly by the plates
1169 theG4AssemblyVolume->AddPlacedVolume(logvol, transl, rotmat);
1170#ifdef G4VERBOSE
1172 {
1173 G4cout << " G4AssemblyVolume->AddPlacedVolume " << ii << " "
1174 << logvol->GetName() << " translation " << transl
1175 << " rotmat " << rotmat->colX() << " " << rotmat->colY()
1176 << " " << rotmat->colZ() << G4endl;
1177 }
1178#endif
1179 }
1180 }
1181
1182 // Rotation and Translation of the assembly inside the world
1183
1184 G4tgrPlaceSimple* placeSimple = (G4tgrPlaceSimple*) place;
1185 G4String rmName = placeSimple->GetRotMatName();
1186 G4RotationMatrix* rotmat =
1188 G4ThreeVector transl = place->GetPlacement();
1189
1190 G4LogicalVolume* parentLV_nonconst =
1191 const_cast<G4LogicalVolume*>(parentLV);
1192 theG4AssemblyVolume->MakeImprint(parentLV_nonconst, transl, rotmat);
1193 }
1194 else // If it is G4tgrVolumeAssembly
1195 {
1196 G4String ErrMessage =
1197 "Volume type not supported: " + theTgrVolume->GetType() + ", sorry...";
1198 G4Exception("G4tgbVolume::ConstructG4PhysVol()", "NotImplemented",
1199 FatalException, ErrMessage);
1200 }
1201 }
1202
1203 return physvol;
1204}
CLHEP::Hep3Vector G4ThreeVector
HepGeom::Rotate3D G4Rotate3D
G4GLOB_DLL 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)
G4LogicalVolume * FindG4LogVol(const G4String &theName, const G4bool bExists=false)
G4tgbVolume * FindVolume(const G4String &volname)
static G4tgbVolumeMgr * GetInstance()
G4VSolid * BuildSolidForDivision(G4VSolid *parentSolid, EAxis axis)
void ConstructG4Volumes(const G4tgrPlace *place, const G4LogicalVolume *parentLV)
Definition: G4tgbVolume.cc:121
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:44
const G4String & GetType() const
Definition: G4tgrPlace.hh:55
unsigned int GetCopyNo() const
Definition: G4tgrPlace.hh:54
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:94
const G4String & GetType() const
Definition: G4tgrVolume.hh:85
EAxis
Definition: geomdefs.hh:54

Referenced by ConstructG4Volumes().

◆ ConstructG4Volumes()

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

Definition at line 121 of file G4tgbVolume.cc.

123{
124#ifdef G4VERBOSE
126 {
127 G4cout << G4endl << "@@@ G4tgbVolume::ConstructG4Volumes - " << GetName()
128 << G4endl;
129 if(place && parentLV)
130 G4cout << " place in LV " << parentLV->GetName() << G4endl;
131 }
132#endif
134 G4LogicalVolume* logvol = g4vmgr->FindG4LogVol(GetName());
135 G4bool bFirstCopy = false;
136 G4VPhysicalVolume* physvol = nullptr;
137 if(logvol == nullptr)
138 {
139 bFirstCopy = true;
140 if(theTgrVolume->GetType() != "VOLDivision")
141 {
142 //--- If first time build solid and LogVol
143 G4VSolid* solid = FindOrConstructG4Solid(theTgrVolume->GetSolid());
144 if(solid != nullptr) // for G4AssemblyVolume it is nullptr
145 {
146 g4vmgr->RegisterMe(solid);
147 logvol = ConstructG4LogVol(solid);
148 g4vmgr->RegisterMe(logvol);
149 g4vmgr->RegisterChildParentLVs(logvol, parentLV);
150 }
151 }
152 else
153 {
154 return;
155 }
156 }
157 //--- Construct PhysVol
158 physvol = ConstructG4PhysVol(place, logvol, parentLV);
159
160 if(physvol != nullptr) // nullptr for G4AssemblyVolumes
161 {
162 g4vmgr->RegisterMe(physvol);
163
164 if(logvol == nullptr) // case of divisions
165 {
166 logvol = physvol->GetLogicalVolume();
167 }
168 }
169 else
170 {
171 return;
172 }
173
174 //--- If first copy build children placements in this LogVol
175 if(bFirstCopy)
176 {
177 std::pair<G4mmapspl::iterator, G4mmapspl::iterator> children =
179 for(auto cite = children.first; cite != children.second; ++cite)
180 {
181 //----- Call G4tgrPlace ->constructG4Volumes
182 //---- find G4tgbVolume corresponding to the G4tgrVolume
183 // pointed by G4tgrPlace
184 G4tgrPlace* pl = const_cast<G4tgrPlace*>((*cite).second);
185 G4tgbVolume* svol = g4vmgr->FindVolume(pl->GetVolume()->GetName());
186 //--- find copyNo
187#ifdef G4VERBOSE
189 {
190 G4cout << " G4tgbVolume::ConstructG4Volumes - construct daughter "
191 << pl->GetVolume()->GetName() << " # " << pl->GetCopyNo()
192 << G4endl;
193 }
194#endif
195 svol->ConstructG4Volumes(pl, logvol);
196 }
197 }
198}
bool G4bool
Definition: G4Types.hh:86
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:854
G4VSolid * FindOrConstructG4Solid(const G4tgrSolid *vol)
Definition: G4tgbVolume.cc:201
G4LogicalVolume * ConstructG4LogVol(const G4VSolid *solid)
Definition: G4tgbVolume.cc:765
G4tgrVolume * GetVolume() const
Definition: G4tgrPlace.hh:53
static G4tgrVolumeMgr * GetInstance()
std::pair< G4mmapspl::iterator, G4mmapspl::iterator > GetChildren(const G4String &name)
const G4String & GetName() const
Definition: G4tgrVolume.hh:83
G4tgrSolid * GetSolid() const
Definition: G4tgrVolume.hh:86

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

◆ FindOrConstructG4Solid()

G4VSolid * G4tgbVolume::FindOrConstructG4Solid ( const G4tgrSolid vol)

Definition at line 201 of file G4tgbVolume.cc.

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

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

Referenced by ConstructG4LogVol().

◆ GetName()

const G4String & G4tgbVolume::GetName ( ) const
inline

Definition at line 99 of file G4tgbVolume.hh.

99{ return theTgrVolume->GetName(); }

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

◆ GetVisibility()

G4bool G4tgbVolume::GetVisibility ( ) const
inline

Definition at line 100 of file G4tgbVolume.hh.

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

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: