Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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::G4tgbVolume ( )

Definition at line 109 of file G4tgbVolume.cc.

110  : theTgrVolume(0), theG4AssemblyVolume(0)
111 {
112 }
G4tgbVolume::~G4tgbVolume ( )

Definition at line 116 of file G4tgbVolume.cc.

117 {
118 }
G4tgbVolume::G4tgbVolume ( G4tgrVolume vol)

Definition at line 122 of file G4tgbVolume.cc.

123 {
124  theTgrVolume = vol;
125  theG4AssemblyVolume = 0;
126 }

Member Function Documentation

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

Definition at line 1183 of file G4tgbVolume.cc.

1184 {
1185  G4VSolid* solid=0;
1186  G4double redf = (parentSolid->GetExtent().GetXmax()-parentSolid->GetExtent().GetXmin());
1187  redf = std::min(redf,parentSolid->GetExtent().GetYmax()-parentSolid->GetExtent().GetYmin());
1188  redf = std::min(redf,parentSolid->GetExtent().GetZmax()-parentSolid->GetExtent().GetZmin());
1189  redf *= 0.001; //make daugther much smaller, to fit in parent
1190 
1191  if( parentSolid->GetEntityType() == "G4Box" )
1192  {
1193  G4Box* psolid = (G4Box*)(parentSolid);
1194  solid = new G4Box(GetName(), psolid->GetXHalfLength()*redf,
1195  psolid->GetZHalfLength()*redf,
1196  psolid->GetZHalfLength()*redf);
1197  }
1198  else if ( parentSolid->GetEntityType() == "G4Tubs" )
1199  {
1200  G4Tubs* psolid = (G4Tubs*)(parentSolid);
1201  solid = new G4Tubs( GetName(), psolid->GetInnerRadius()*redf,
1202  psolid->GetOuterRadius()*redf,
1203  psolid->GetZHalfLength()*redf,
1204  psolid->GetStartPhiAngle(),
1205  psolid->GetDeltaPhiAngle());
1206  }
1207  else if ( parentSolid->GetEntityType() == "G4Cons" )
1208  {
1209  G4Cons* psolid = (G4Cons*)(parentSolid);
1210  solid = new G4Cons( GetName(), psolid->GetInnerRadiusMinusZ()*redf,
1211  psolid->GetOuterRadiusMinusZ()*redf,
1212  psolid->GetInnerRadiusPlusZ()*redf,
1213  psolid->GetOuterRadiusPlusZ()*redf,
1214  psolid->GetZHalfLength()*redf,
1215  psolid->GetStartPhiAngle(),
1216  psolid->GetDeltaPhiAngle());
1217  }
1218  else if ( parentSolid->GetEntityType() == "G4Trd" )
1219  {
1220  G4Trd* psolid = (G4Trd*)(parentSolid);
1221  G4double mpDx1 = psolid->GetXHalfLength1();
1222  G4double mpDx2 = psolid->GetXHalfLength2();
1223 
1224  if( axis == kXAxis && std::fabs(mpDx1 - mpDx2) > G4GeometryTolerance::GetInstance()->GetSurfaceTolerance() )
1225  {
1226  solid = new G4Trap( GetName(), psolid->GetZHalfLength()*redf,
1227  psolid->GetYHalfLength1()*redf,
1228  psolid->GetXHalfLength2()*redf,
1229  psolid->GetXHalfLength1()*redf );
1230  }
1231  else
1232  {
1233  solid = new G4Trd( GetName(), psolid->GetXHalfLength1()*redf,
1234  psolid->GetXHalfLength2()*redf,
1235  psolid->GetYHalfLength1()*redf,
1236  psolid->GetYHalfLength2()*redf,
1237  psolid->GetZHalfLength()*redf);
1238  }
1239 
1240  }
1241  else if ( parentSolid->GetEntityType() == "G4Para" )
1242  {
1243  G4Para* psolid = (G4Para*)(parentSolid);
1244  solid = new G4Para( GetName(), psolid->GetXHalfLength()*redf,
1245  psolid->GetYHalfLength()*redf,
1246  psolid->GetZHalfLength()*redf,
1247  std::atan(psolid->GetTanAlpha()),
1248  psolid->GetSymAxis().theta(),
1249  psolid->GetSymAxis().phi() );
1250  }
1251  else if ( parentSolid->GetEntityType() == "G4Polycone" )
1252  {
1253  G4Polycone* psolid = (G4Polycone*)(parentSolid);
1254  G4PolyconeHistorical origParam = *(psolid->GetOriginalParameters());
1255  for( G4int ii = 0; ii < origParam.Num_z_planes; ii++ )
1256  {
1257  origParam.Rmin[ii] = origParam.Rmin[ii]*redf;
1258  origParam.Rmax[ii] = origParam.Rmax[ii]*redf;
1259  }
1260  solid = new G4Polycone( GetName(), psolid->GetStartPhi(),
1261  psolid->GetEndPhi(),
1262  origParam.Num_z_planes, origParam.Z_values,
1263  origParam.Rmin, origParam.Rmax);
1264 
1265  }
1266  else if ( parentSolid->GetEntityType() == "G4GenericPolycone" )
1267  {
1268  G4GenericPolycone* psolid = (G4GenericPolycone*)(parentSolid);
1269  const G4int numRZ = psolid->GetNumRZCorner();
1270  G4double* r = new G4double[numRZ];
1271  G4double* z = new G4double[numRZ];
1272  for( G4int ii = 0; ii < numRZ; ii++ )
1273  {
1274  r[ii] = psolid->GetCorner(ii).r;
1275  z[ii] = psolid->GetCorner(ii).z;
1276  }
1277  solid = new G4GenericPolycone( GetName(), psolid->GetStartPhi(),
1278  psolid->GetEndPhi() - psolid->GetStartPhi(),
1279  numRZ, r, z);
1280  delete [] r; delete [] z;
1281 
1282  }
1283  else if ( parentSolid->GetEntityType() == "G4Polyhedra" )
1284  {
1285  G4Polyhedra* psolid = (G4Polyhedra*)(parentSolid);
1286  G4PolyhedraHistorical origParam = *(psolid->GetOriginalParameters());
1287  for( G4int ii = 0; ii < origParam.Num_z_planes; ii++ )
1288  {
1289  origParam.Rmin[ii] = origParam.Rmin[ii]*redf;
1290  origParam.Rmax[ii] = origParam.Rmax[ii]*redf;
1291  }
1292  solid = new G4Polyhedra( GetName(), psolid->GetStartPhi(),
1293  psolid->GetEndPhi(),
1294  psolid->GetNumSide(),
1295  origParam.Num_z_planes, origParam.Z_values,
1296  origParam.Rmin, origParam.Rmax);
1297  }
1298  else
1299  {
1300  G4String ErrMessage = "Solid type not supported. VOLUME= " + GetName()
1301  + " Solid type= " + parentSolid->GetEntityType() + "\n"
1302  + "Only supported types are: G4Box, G4Tubs, G4Cons,"
1303  + " G4Trd, G4Para, G4Polycone, G4Polyhedra.";
1304  G4Exception("G4tgbVolume::BuildSolidForDivision()", "NotImplemented",
1305  FatalException, ErrMessage);
1306  return 0;
1307  }
1308 
1309 #ifdef G4VERBOSE
1310  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1311  {
1312  G4cout << " Constructing new G4Solid for division: "
1313  << *solid << G4endl;
1314  }
1315 #endif
1316  return solid;
1317 }
G4ThreeVector GetSymAxis() const
G4PolyconeSideRZ GetCorner(G4int index) const
G4double GetXHalfLength() const
Definition: G4Para.hh:77
G4double GetXmin() const
Definition: G4VisExtent.hh:89
G4double GetYHalfLength1() const
Definition: G4Box.hh:64
Definition: G4Tubs.hh:85
G4double GetXmax() const
Definition: G4VisExtent.hh:90
const G4String & GetName() const
Definition: G4tgbVolume.hh:108
G4double GetOuterRadiusMinusZ() const
Definition: G4Trd.hh:72
virtual G4GeometryType GetEntityType() const =0
G4double GetZHalfLength() const
G4double GetEndPhi() const
G4double GetEndPhi() const
int G4int
Definition: G4Types.hh:78
G4double GetZHalfLength() const
G4double GetEndPhi() const
G4double GetXHalfLength2() const
G4int GetNumSide() const
G4GLOB_DLL std::ostream G4cout
G4double GetDeltaPhiAngle() const
static G4int GetVerboseLevel()
G4double GetYmax() const
Definition: G4VisExtent.hh:92
Definition: G4Cons.hh:83
G4double GetStartPhiAngle() const
G4double GetYHalfLength2() const
G4double GetZmax() const
Definition: G4VisExtent.hh:94
G4double GetStartPhiAngle() const
double phi() const
G4double GetStartPhi() const
G4double GetInnerRadiusPlusZ() const
G4double GetInnerRadius() const
G4PolyconeHistorical * GetOriginalParameters() const
G4double GetStartPhi() const
G4double GetXHalfLength() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double theta() const
virtual G4VisExtent GetExtent() const
Definition: G4VSolid.cc:642
G4double GetTanAlpha() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4PolyhedraHistorical * GetOriginalParameters() const
G4double GetZmin() const
Definition: G4VisExtent.hh:93
G4double GetZHalfLength() const
G4double GetZHalfLength() const
tuple z
Definition: test.py:28
G4double GetYmin() const
Definition: G4VisExtent.hh:91
G4double GetStartPhi() const
#define G4endl
Definition: G4ios.hh:61
G4double GetXHalfLength1() const
double G4double
Definition: G4Types.hh:76
G4double GetInnerRadiusMinusZ() const
G4int GetNumRZCorner() const
G4double GetOuterRadiusPlusZ() const
G4double GetZHalfLength() const
static G4GeometryTolerance * GetInstance()
G4double GetOuterRadius() const
G4double GetDeltaPhiAngle() const
G4double GetYHalfLength() const

Here is the caller graph for this function:

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

Definition at line 705 of file G4tgbVolume.cc.

708 {
709  if( NoParamExpected != NoParam )
710  {
711  G4String Err1 = "Solid type " + solidType + " should have ";
712  G4String Err2 = G4UIcommand::ConvertToString(G4int(NoParamExpected))
713  + " parameters,\n";
714  G4String Err3 = "and it has "
716  G4String ErrMessage = Err1 + Err2 + Err3 + " !";
717  G4Exception("G4tgbVolume::CheckNoSolidParams()", "InvalidSetup",
718  FatalException, ErrMessage);
719  }
720 }
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:372
int G4int
Definition: G4Types.hh:78
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

Here is the caller graph for this function:

G4LogicalVolume * G4tgbVolume::ConstructG4LogVol ( const G4VSolid solid)

Definition at line 724 of file G4tgbVolume.cc.

725 {
726  G4LogicalVolume* logvol;
727 
728 #ifdef G4VERBOSE
730  {
731  G4cout << " G4tgbVolume::ConstructG4LogVol() - " << GetName() << G4endl;
732  }
733 #endif
734 
735  //----------- Get the material first
737  ->FindOrBuildG4Material( theTgrVolume->GetMaterialName() );
738  if( mate == 0 )
739  {
740  G4String ErrMessage = "Material not found "
741  + theTgrVolume->GetMaterialName()
742  + " for volume " + GetName() + ".";
743  G4Exception("G4tgbVolume::ConstructG4LogVol()", "InvalidSetup",
744  FatalException, ErrMessage);
745  }
746 #ifdef G4VERBOSE
748  {
749  G4cout << " G4tgbVolume::ConstructG4LogVol() -"
750  << " Material constructed: " << mate->GetName() << G4endl;
751  }
752 #endif
753 
754  //---------- Construct the LV
755  logvol = new G4LogicalVolume( const_cast<G4VSolid*>(solid),
756  const_cast<G4Material*>(mate), GetName() );
757 
758 #ifdef G4VERBOSE
760  {
761  G4cout << " Constructing new G4LogicalVolume: "
762  << logvol->GetName() << " mate " << mate->GetName() << G4endl;
763  }
764 #endif
765 
766  //---------- Set visibility and colour
767  if( !GetVisibility() || GetColour()[0] != -1 )
768  {
769  G4VisAttributes* visAtt = new G4VisAttributes();
770 #ifdef G4VERBOSE
772  {
773  G4cout << " Constructing new G4VisAttributes: "
774  << *visAtt << G4endl;
775  }
776 #endif
777 
778  if( !GetVisibility() )
779  {
780  visAtt->SetVisibility( GetVisibility() );
781  }
782  else if( GetColour()[0] != -1 )
783  {
784  // this else should not be necessary, because if the visibility
785  // is set to off, colour should have no effect. But it does not
786  // work: if you set colour and vis off, it is visualized!?!?!?
787 
788  const G4double* col = GetColour();
789  if( col[3] == -1. )
790  {
791  visAtt->SetColour( G4Colour(col[0],col[1],col[2]));
792  }
793  else
794  {
795  visAtt->SetColour( G4Colour(col[0],col[1],col[2],col[3]));
796  }
797  }
798  logvol->SetVisAttributes(visAtt);
799  }
800 
801 #ifdef G4VERBOSE
803  {
804  G4cout << " G4tgbVolume::ConstructG4LogVol() -"
805  << " Created logical volume: " << GetName() << G4endl;
806  }
807 #endif
808 
809  return logvol;
810 }
void SetColour(const G4Colour &)
const G4String & GetName() const
Definition: G4Material.hh:178
G4bool GetVisibility() const
Definition: G4tgbVolume.hh:109
const G4double * GetColour() const
Definition: G4tgbVolume.hh:110
const G4String & GetName() const
Definition: G4tgbVolume.hh:108
G4GLOB_DLL std::ostream G4cout
static G4int GetVerboseLevel()
G4Material * FindOrBuildG4Material(const G4String &name, G4bool bMustExist=1)
void SetVisibility(G4bool=true)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
const G4String & GetName() const
double G4double
Definition: G4Types.hh:76
static G4tgbMaterialMgr * GetInstance()
const G4String & GetMaterialName() const
Definition: G4tgrVolume.hh:93
void SetVisAttributes(const G4VisAttributes *pVA)

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 815 of file G4tgbVolume.cc.

818 {
819  G4VPhysicalVolume* physvol = 0;
820  G4int copyNo;
821 
822  //----- Case of placement of top volume
823  if( place == 0 )
824  {
825 #ifdef G4VERBOSE
827  {
828  G4cout << " G4tgbVolume::ConstructG4PhysVol() - World: "
829  << GetName() << G4endl;
830  }
831 #endif
832  physvol = new G4PVPlacement(0, G4ThreeVector(),
833  const_cast<G4LogicalVolume*>(currentLV),
834  GetName(), 0, false, 0,
835  theTgrVolume->GetCheckOverlaps());
836 #ifdef G4VERBOSE
838  {
839  G4cout << " Constructing new : G4PVPlacement "
840  << physvol->GetName() << G4endl;
841  }
842 #endif
843  }
844  else
845  {
846  copyNo = place->GetCopyNo();
847 
848 #ifdef G4VERBOSE
850  {
851  G4cout << " G4tgbVolume::ConstructG4PhysVol() - " << GetName() << G4endl
852  << " inside " << parentLV->GetName() << " copy No: " << copyNo
853  << " of type= " << theTgrVolume->GetType() << G4endl
854  << " placement type= " << place->GetType() << G4endl;
855  }
856 #endif
857 
858  if( theTgrVolume->GetType() == "VOLSimple" )
859  {
860  //----- Get placement
861 #ifdef G4VERBOSE
863  {
864  G4cout << " G4tgbVolume::ConstructG4PhysVol() - Placement type = "
865  << place->GetType() << G4endl;
866  }
867 #endif
868 
869  //--------------- If it is G4tgrPlaceSimple
870  if( place->GetType() == "PlaceSimple" )
871  {
872  //----- Get rotation matrix
873  G4tgrPlaceSimple* placeSimple = (G4tgrPlaceSimple*)place;
874  G4String rmName = placeSimple->GetRotMatName();
875 
877  ->FindOrBuildG4RotMatrix( rmName );
878  //----- Place volume in mother
879  G4double check = (rotmat->colX().cross(rotmat->colY()))*rotmat->colZ();
880  G4double tol = 1.0e-3;
881  //---- Check that matrix is ortogonal
882  if (1-std::abs(check)>tol)
883  {
884  G4cerr << " Matrix : " << rmName << " " << rotmat->colX()
885  << " " << rotmat->colY() << " " << rotmat->colZ() << G4endl
886  << " product x X y * z = " << check << " x X y "
887  << rotmat->colX().cross(rotmat->colY()) << G4endl;
888  G4String ErrMessage = "Rotation is not ortogonal " + rmName + " !";
889  G4Exception("G4tgbVolume::ConstructG4PhysVol()",
890  "InvalidSetup", FatalException, ErrMessage);
891  //---- Check if it is reflection
892  }
893  else if (1+check<=tol)
894  {
895  G4Translate3D transl = place->GetPlacement();
896  G4Transform3D trfrm = transl * G4Rotate3D(*rotmat);
897  physvol = (G4ReflectionFactory::Instance()->Place(trfrm, GetName(),
898  const_cast<G4LogicalVolume*>(currentLV),
899  const_cast<G4LogicalVolume*>(parentLV),
900  false, copyNo, false )).first;
901  }
902  else
903  {
904 #ifdef G4VERBOSE
906  {
907  G4cout << "Construction new G4VPhysicalVolume"
908  << " through G4ReflectionFactory " << GetName()
909  << " in volume " << parentLV->GetName()
910  << " copyNo " << copyNo
911  << " position " << place->GetPlacement()
912  << " ROT " << rotmat->colX()
913  << " " << rotmat->colY()
914  << " " << rotmat->colZ() << G4endl;
915  }
916 #endif
917  physvol = new G4PVPlacement( rotmat, place->GetPlacement(),
918  const_cast<G4LogicalVolume*>(currentLV),
919  GetName(),
920  const_cast<G4LogicalVolume*>(parentLV),
921  false, copyNo,
922  theTgrVolume->GetCheckOverlaps());
923  }
924 
925  //--------------- If it is G4tgrPlaceParam
926  }
927  else if( place->GetType() == "PlaceParam" )
928  {
930 
931  //----- See what parameterisation type
932 #ifdef G4VERBOSE
934  {
935  G4cout << " G4tgbVolume::ConstructG4PhysVol() -" << G4endl
936  << " param: " << GetName() << " in " << parentLV->GetName()
937  << " param type= " << dp->GetParamType() << G4endl;
938  }
939 #endif
940 
941  G4tgbPlaceParameterisation * param=0;
942 
943  if( (dp->GetParamType() == "CIRCLE")
944  || (dp->GetParamType() == "CIRCLE_XY")
945  || (dp->GetParamType() == "CIRCLE_XZ")
946  || (dp->GetParamType() == "CIRCLE_YZ") )
947  {
948  param = new G4tgbPlaceParamCircle(dp);
949 
950  }
951  else if( (dp->GetParamType() == "LINEAR")
952  || (dp->GetParamType() == "LINEAR_X")
953  || (dp->GetParamType() == "LINEAR_Y")
954  || (dp->GetParamType() == "LINEAR_Z") )
955  {
956  param = new G4tgbPlaceParamLinear(dp);
957 
958  }
959  else if( (dp->GetParamType() == "SQUARE")
960  || (dp->GetParamType() == "SQUARE_XY")
961  || (dp->GetParamType() == "SQUARE_XZ")
962  || (dp->GetParamType() == "SQUARE_YZ") )
963  {
964  param = new G4tgbPlaceParamSquare(dp);
965  }
966  else
967  {
968  G4String ErrMessage = "Parameterisation has wrong type, TYPE: "
969  + G4String(dp->GetParamType()) + " !";
970  G4Exception("G4tgbVolume::ConstructG4PhysVol", "WrongArgument",
971  FatalException, ErrMessage);
972  return 0;
973  }
974 #ifdef G4VERBOSE
976  {
977  G4cout << " G4tgbVolume::ConstructG4PhysVol() -" << G4endl
978  << " New G4PVParameterised: " << GetName() << " vol "
979  << currentLV->GetName() << " in vol " << parentLV->GetName()
980  << " axis " << param->GetAxis() << " nCopies "
981  << param->GetNCopies() << G4endl;
982  }
983 #endif
984  physvol = new G4PVParameterised(GetName(),
985  const_cast<G4LogicalVolume*>(currentLV),
986  const_cast<G4LogicalVolume*>(parentLV),
987  EAxis(param->GetAxis()),
988  param->GetNCopies(), param);
989 #ifdef G4VERBOSE
991  {
992  G4cout << " Constructing new G4PVParameterised: "
993  << physvol->GetName() << " in volume " << parentLV->GetName()
994  << " N copies " << param->GetNCopies()
995  << " axis " << param->GetAxis() << G4endl;
996  }
997 #endif
998 
999  }
1000  else if( place->GetType() == "PlaceReplica" )
1001  {
1002  //--------------- If it is PlaceReplica
1003  G4tgrPlaceDivRep* dpr = (G4tgrPlaceDivRep*)place;
1004 
1005 #ifdef G4VERBOSE
1006  if( G4tgrMessenger::GetVerboseLevel() >= 2 )
1007  {
1008  G4cout << " G4tgbVolume::ConstructG4PhysVol() -" << G4endl
1009  << " replica" << " " << currentLV->GetName()
1010  << " in " << parentLV->GetName()
1011  << " NDiv " << dpr->GetNDiv() << " Width " << dpr->GetWidth()
1012  << " offset " << dpr->GetOffset() << G4endl;
1013  }
1014 #endif
1015  physvol = new G4PVReplica(GetName(),
1016  const_cast<G4LogicalVolume*>(currentLV),
1017  const_cast<G4LogicalVolume*>(parentLV),
1018  EAxis(dpr->GetAxis()), dpr->GetNDiv(),
1019  dpr->GetWidth(), dpr->GetOffset());
1020 #ifdef G4VERBOSE
1021  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1022  {
1023  G4cout << " Constructing new G4PVReplica: "
1024  << currentLV->GetName()
1025  << " in " << parentLV->GetName()
1026  << " NDiv " << dpr->GetNDiv() << " Width " << dpr->GetWidth()
1027  << " offset " << dpr->GetOffset() << G4endl;
1028  }
1029 #endif
1030  }
1031  }
1032  else if( theTgrVolume->GetType() == "VOLDivision" )
1033  {
1034  G4tgrVolumeDivision* volr = (G4tgrVolumeDivision*)theTgrVolume;
1035  G4tgrPlaceDivRep* placeDiv = volr->GetPlaceDivision() ;
1036  G4VSolid* solid = BuildSolidForDivision( parentLV->GetSolid(), placeDiv->GetAxis() );
1038  ->FindOrBuildG4Material( theTgrVolume->GetMaterialName() );
1039  G4LogicalVolume* divLV = new G4LogicalVolume(solid,
1040  const_cast<G4Material*>(mate),
1041  GetName() );
1042 #ifdef G4VERBOSE
1043  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1044  {
1045  G4cout << " Constructed new G4LogicalVolume for division: "
1046  << divLV->GetName() << " mate " << mate->GetName() << G4endl;
1047  }
1048 #endif
1049 
1050  G4DivType divType = placeDiv->GetDivType();
1051  switch (divType)
1052  {
1053  case DivByNdiv:
1054  physvol = new G4PVDivision(GetName(), (G4LogicalVolume*)divLV,
1055  const_cast<G4LogicalVolume*>(parentLV),
1056  placeDiv->GetAxis(), placeDiv->GetNDiv(),
1057  placeDiv->GetOffset());
1058 #ifdef G4VERBOSE
1059  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1060  {
1061  G4cout << " Constructing new G4PVDivision by number of divisions: "
1062  << GetName() << " in " << parentLV->GetName()
1063  << " axis " << placeDiv->GetAxis()
1064  << " Ndiv " << placeDiv->GetNDiv()
1065  << " offset " << placeDiv->GetOffset() << G4endl;
1066  }
1067 #endif
1068  break;
1069  case DivByWidth:
1070  physvol = new G4PVDivision(GetName(), (G4LogicalVolume*)divLV,
1071  const_cast<G4LogicalVolume*>(parentLV),
1072  placeDiv->GetAxis(), placeDiv->GetWidth(),
1073  placeDiv->GetOffset());
1074 #ifdef G4VERBOSE
1075  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1076  {
1077  G4cout << " Constructing new G4PVDivision by width: "
1078  << GetName() << " in " << parentLV->GetName()
1079  << " axis " << placeDiv->GetAxis()
1080  << " width " << placeDiv->GetWidth()
1081  << " offset " << placeDiv->GetOffset() << G4endl;
1082  }
1083 #endif
1084  break;
1085  case DivByNdivAndWidth:
1086  physvol = new G4PVDivision(GetName(), (G4LogicalVolume*)divLV,
1087  const_cast<G4LogicalVolume*>(parentLV),
1088  placeDiv->GetAxis(), placeDiv->GetNDiv(),
1089  placeDiv->GetWidth(),
1090  placeDiv->GetOffset());
1091 #ifdef G4VERBOSE
1092  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1093  {
1094  G4cout << " Constructing new G4PVDivision by width"
1095  << " and number of divisions: "
1096  << GetName() << " in " << parentLV->GetName()
1097  << " axis " << placeDiv->GetAxis()
1098  << " Ndiv " << placeDiv->GetNDiv()
1099  << " width " << placeDiv->GetWidth()
1100  << " offset " << placeDiv->GetOffset() << G4endl;
1101  }
1102 #endif
1103  break;
1104  }
1105  }
1106  else if( theTgrVolume->GetType() == "VOLAssembly" )
1107  {
1108  // Define one layer as one assembly volume
1109  G4tgrVolumeAssembly * tgrAssembly = (G4tgrVolumeAssembly *)theTgrVolume;
1110 
1111  if( !theG4AssemblyVolume )
1112  {
1113  theG4AssemblyVolume = new G4AssemblyVolume;
1114 #ifdef G4VERBOSE
1115  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1116  {
1117  G4cout << " Constructing new G4AssemblyVolume: "
1118  << " number of assembly components "
1119  << tgrAssembly->GetNoComponents() << G4endl;
1120  }
1121 #endif
1123  for( G4int ii = 0; ii < tgrAssembly->GetNoComponents(); ii++ )
1124  {
1125  // Rotation and translation of a plate inside the assembly
1126 
1127  G4ThreeVector transl = tgrAssembly->GetComponentPos(ii);
1128  G4String rmName = tgrAssembly->GetComponentRM(ii);
1130  ->FindOrBuildG4RotMatrix( rmName );
1131 
1132  //----- Get G4LogicalVolume of component
1133  G4String lvname = tgrAssembly->GetComponentName(ii);
1134  G4LogicalVolume* logvol = g4vmgr->FindG4LogVol( lvname);
1135  if( logvol == 0 )
1136  {
1137  g4vmgr->FindVolume( lvname )->ConstructG4Volumes( 0, 0);
1138  logvol = g4vmgr->FindG4LogVol( lvname, true );
1139  }
1140  // Fill the assembly by the plates
1141  theG4AssemblyVolume->AddPlacedVolume( logvol, transl, rotmat );
1142 #ifdef G4VERBOSE
1143  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1144  {
1145  G4cout << " G4AssemblyVolume->AddPlacedVolume " << ii
1146  << " " << logvol->GetName()
1147  << " translation " << transl
1148  << " rotmat " << rotmat->colX()
1149  << " " << rotmat->colY()
1150  << " " << rotmat->colZ() << G4endl;
1151  }
1152 #endif
1153  }
1154  }
1155 
1156  // Rotation and Translation of the assembly inside the world
1157 
1158  G4tgrPlaceSimple* placeSimple = (G4tgrPlaceSimple*)place;
1159  G4String rmName = placeSimple->GetRotMatName();
1161  ->FindOrBuildG4RotMatrix( rmName );
1162  G4ThreeVector transl = place->GetPlacement();
1163 
1164  G4LogicalVolume* parentLV_nonconst =
1165  const_cast<G4LogicalVolume*>(parentLV);
1166  theG4AssemblyVolume->MakeImprint( parentLV_nonconst, transl, rotmat );
1167 
1168  }
1169  else // If it is G4tgrVolumeAssembly
1170  {
1171  G4String ErrMessage = "Volume type not supported: "
1172  + theTgrVolume->GetType() + ", sorry...";
1173  G4Exception("G4tgbVolume::ConstructG4PhysVol()", "NotImplemented",
1174  FatalException, ErrMessage);
1175  }
1176  }
1177 
1178  return physvol;
1179 }
const G4String & GetComponentName(G4int ii) const
G4bool GetCheckOverlaps() const
Definition: G4tgrVolume.hh:100
void MakeImprint(G4LogicalVolume *pMotherLV, G4ThreeVector &translationInMother, G4RotationMatrix *pRotationInMother, G4int copyNumBase=0, G4bool surfCheck=false)
CLHEP::Hep3Vector G4ThreeVector
Hep3Vector colX() const
G4int GetNoComponents() const
static G4tgbRotationMatrixMgr * GetInstance()
const G4String & GetType() const
Definition: G4tgrPlace.hh:61
G4VSolid * GetSolid() const
G4tgrPlaceDivRep * GetPlaceDivision()
const G4String & GetName() const
Definition: G4tgbVolume.hh:108
int G4int
Definition: G4Types.hh:78
static G4ReflectionFactory * Instance()
virtual G4ThreeVector GetPlacement() const
Definition: G4tgrPlace.cc:52
G4GLOB_DLL std::ostream G4cout
G4PhysicalVolumesPair Place(const G4Transform3D &transform3D, const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, G4bool isMany, G4int copyNo, G4bool surfCheck=false)
const G4String & GetName() const
static G4int GetVerboseLevel()
G4ThreeVector GetComponentPos(G4int ii) const
G4Material * FindOrBuildG4Material(const G4String &name, G4bool bMustExist=1)
G4VSolid * BuildSolidForDivision(G4VSolid *parentSolid, EAxis axis)
EAxis GetAxis() const
G4DivType
const G4String & GetRotMatName() const
Hep3Vector colZ() const
const G4String & GetParamType() const
Hep3Vector colY() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4DivType GetDivType() const
const G4String & GetType() const
Definition: G4tgrVolume.hh:91
HepGeom::Rotate3D G4Rotate3D
const G4String & GetComponentRM(G4int ii) const
void AddPlacedVolume(G4LogicalVolume *pPlacedVolume, G4ThreeVector &translation, G4RotationMatrix *rotation)
G4double GetOffset() const
EAxis
Definition: geomdefs.hh:54
unsigned int GetCopyNo() const
Definition: G4tgrPlace.hh:60
G4tgbVolume * FindVolume(const G4String &volname)
void ConstructG4Volumes(const G4tgrPlace *place, const G4LogicalVolume *parentLV)
Definition: G4tgbVolume.cc:130
G4RotationMatrix * FindOrBuildG4RotMatrix(const G4String &name)
#define G4endl
Definition: G4ios.hh:61
Hep3Vector cross(const Hep3Vector &) const
const G4String & GetName() const
double G4double
Definition: G4Types.hh:76
static G4tgbVolumeMgr * GetInstance()
static G4tgbMaterialMgr * GetInstance()
const G4String & GetMaterialName() const
Definition: G4tgrVolume.hh:93
G4double GetWidth() const
G4GLOB_DLL std::ostream G4cerr
G4int GetNDiv() const
G4LogicalVolume * FindG4LogVol(const G4String &theName, const G4bool bExists=0)

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 130 of file G4tgbVolume.cc.

132 {
133 #ifdef G4VERBOSE
135  {
136  G4cout << G4endl << "@@@ G4tgbVolume::ConstructG4Volumes - " << GetName() << G4endl;
137  if( place && parentLV ) G4cout << " place in LV " << parentLV->GetName() << G4endl;
138  }
139 #endif
141  G4LogicalVolume* logvol = g4vmgr->FindG4LogVol( GetName() );
142  G4bool bFirstCopy = false;
143  if( logvol == 0 )
144  {
145  bFirstCopy = true;
146  if( theTgrVolume->GetType() != "VOLDivision" )
147  {
148  //--- If first time build solid and LogVol
149  G4VSolid* solid = FindOrConstructG4Solid( theTgrVolume->GetSolid() );
150  if( solid != 0 ) // for G4AssemblyVolume it is 0
151  {
152  g4vmgr->RegisterMe( solid );
153  logvol = ConstructG4LogVol( solid );
154  g4vmgr->RegisterMe( logvol );
155  g4vmgr->RegisterChildParentLVs( logvol, parentLV );
156  }
157  }
158  else
159  {
160  return;
161  }
162  }
163  //--- Construct PhysVol
164  G4VPhysicalVolume* physvol = ConstructG4PhysVol( place, logvol, parentLV );
165  if( physvol != 0 ) // 0 for G4AssemblyVolumes
166  {
167  g4vmgr->RegisterMe( physvol );
168 
169  if( logvol == 0 ) // case of divisions
170  {
171  logvol = physvol->GetLogicalVolume();
172  }
173  }
174  else
175  {
176  return;
177  }
178 
179  //--- If first copy build children placements in this LogVol
180  if(bFirstCopy)
181  {
182  std::pair<G4mmapspl::iterator, G4mmapspl::iterator> children
184  G4mmapspl::iterator cite;
185  for( cite = children.first; cite != children.second; cite++ )
186  {
187  //----- Call G4tgrPlace ->constructG4Volumes
188  //---- find G4tgbVolume corresponding to the G4tgrVolume
189  // pointed by G4tgrPlace
190  G4tgrPlace* pl = const_cast<G4tgrPlace*>((*cite).second);
191  G4tgbVolume* svol = g4vmgr->FindVolume( pl->GetVolume()->GetName() );
192  //--- find copyNo
193 #ifdef G4VERBOSE
195  {
196  G4cout << " G4tgbVolume::ConstructG4Volumes - construct daughter " << pl->GetVolume()->GetName() << " # " << pl->GetCopyNo() << G4endl;
197  }
198 #endif
199  svol->ConstructG4Volumes( pl, logvol );
200  }
201  }
202 
203 }
G4VSolid * FindOrConstructG4Solid(const G4tgrSolid *vol)
Definition: G4tgbVolume.cc:208
G4tgrVolume * GetVolume() const
Definition: G4tgrPlace.hh:59
G4LogicalVolume * ConstructG4LogVol(const G4VSolid *solid)
Definition: G4tgbVolume.cc:724
const G4String & GetName() const
Definition: G4tgbVolume.hh:108
tuple pl
Definition: readPY.py:5
G4tgrSolid * GetSolid() const
Definition: G4tgrVolume.hh:92
G4GLOB_DLL std::ostream G4cout
static G4int GetVerboseLevel()
bool G4bool
Definition: G4Types.hh:79
static G4tgrVolumeMgr * GetInstance()
const G4String & GetType() const
Definition: G4tgrVolume.hh:91
G4VPhysicalVolume * ConstructG4PhysVol(const G4tgrPlace *place, const G4LogicalVolume *currentLV, const G4LogicalVolume *parentLV)
Definition: G4tgbVolume.cc:815
void RegisterMe(const G4tgbVolume *vol)
G4LogicalVolume * GetLogicalVolume() const
void RegisterChildParentLVs(const G4LogicalVolume *logvol, const G4LogicalVolume *parentLV)
unsigned int GetCopyNo() const
Definition: G4tgrPlace.hh:60
G4tgbVolume * FindVolume(const G4String &volname)
void ConstructG4Volumes(const G4tgrPlace *place, const G4LogicalVolume *parentLV)
Definition: G4tgbVolume.cc:130
#define G4endl
Definition: G4ios.hh:61
const G4String & GetName() const
static G4tgbVolumeMgr * GetInstance()
const G4String & GetName() const
Definition: G4tgrVolume.hh:89
std::pair< G4mmapspl::iterator, G4mmapspl::iterator > GetChildren(const G4String &name)
G4LogicalVolume * FindG4LogVol(const G4String &theName, const G4bool bExists=0)

Here is the call graph for this function:

Here is the caller graph for this function:

G4VSolid * G4tgbVolume::FindOrConstructG4Solid ( const G4tgrSolid vol)

Definition at line 208 of file G4tgbVolume.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Here is the call graph for this function:

Here is the caller graph for this function:

const G4String& G4tgbVolume::GetName ( void  ) const
inline

Definition at line 108 of file G4tgbVolume.hh.

108 { return theTgrVolume->GetName(); }
const G4String & GetName() const
Definition: G4tgrVolume.hh:89

Here is the call graph for this function:

Here is the caller graph for this function:

G4bool G4tgbVolume::GetVisibility ( ) const
inline

Definition at line 109 of file G4tgbVolume.hh.

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

Here is the call graph for this function:

Here is the caller graph for this function:

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

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