00001
00002
00003
00004
00005
00006 #include "BDSGlobalConstants.hh"
00007 #include "BDSSpoiler.hh"
00008 #include "G4VisAttributes.hh"
00009 #include "G4LogicalVolume.hh"
00010 #include "G4VPhysicalVolume.hh"
00011 #include "G4PVPlacement.hh"
00012 #include "G4UserLimits.hh"
00013 #include "G4TransportationManager.hh"
00014
00015 #include "G4SDManager.hh"
00016
00017 #include <map>
00018
00019
00020
00021
00022
00023 typedef std::map<G4String,int> LogVolCountMap;
00024 extern LogVolCountMap* LogVolCount;
00025
00026 typedef std::map<G4String,G4LogicalVolume*> LogVolMap;
00027 extern LogVolMap* LogVol;
00028
00029 extern BDSMaterials* theMaterials;
00030
00031
00032 BDSSpoiler::BDSSpoiler (G4String& aName,G4double aLength,G4double bpRad,
00033 G4double xAper,G4double yAper,
00034 G4Material* SpoilerMaterial):
00035 BDSAcceleratorComponent(aName,
00036 aLength,bpRad,xAper,yAper,
00037 SetVisAttributes()),
00038 itsPhysiComp(NULL), itsPhysiComp2(NULL), itsSolidLogVol(NULL),
00039 itsInnerLogVol(NULL), itsVisAttributes(NULL), itsEqRhs(NULL),
00040 itsSpoilerMaterial(SpoilerMaterial)
00041 {
00042
00043 if ( (*LogVolCount)[itsName]==0)
00044 {
00045 itsMarkerLogicalVolume=
00046 new G4LogicalVolume(
00047 new G4Box(itsName,
00048 BDSGlobalConstants::Instance()->GetComponentBoxSize()/2,
00049 BDSGlobalConstants::Instance()->GetComponentBoxSize()/2,
00050 itsLength/2),
00051 theMaterials->GetMaterial(BDSGlobalConstants::Instance()->GetVacuumMaterial()),
00052 itsName);
00053 BuildInnerSpoiler();
00054
00055 (*LogVolCount)[itsName]=1;
00056 (*LogVol)[itsName]=itsMarkerLogicalVolume;
00057 }
00058 else
00059 {
00060 (*LogVolCount)[itsName]++;
00061 itsMarkerLogicalVolume=(*LogVol)[itsName];
00062 }
00063 }
00064
00065
00066 G4VisAttributes* BDSSpoiler::SetVisAttributes()
00067 {
00068 itsVisAttributes=new G4VisAttributes(G4Colour(0.3,0.4,0.2));
00069 return itsVisAttributes;
00070 }
00071
00072
00073 void BDSSpoiler::BuildInnerSpoiler()
00074 {
00075 itsSolidLogVol=
00076 new G4LogicalVolume(new G4Box(itsName+"_solid",
00077 BDSGlobalConstants::Instance()->GetComponentBoxSize()/2,
00078 BDSGlobalConstants::Instance()->GetComponentBoxSize()/2,
00079 itsLength/2),
00080 itsSpoilerMaterial,
00081 itsName+"_solid");
00082
00083 itsInnerLogVol=
00084 new G4LogicalVolume(new G4Box(itsName+"_inner",
00085 itsXAper,
00086 itsYAper,
00087 itsLength/2),
00088 theMaterials->GetMaterial(BDSGlobalConstants::Instance()->GetVacuumMaterial()),
00089 itsName+"_inner");
00090
00091 itsPhysiComp2 =
00092 new G4PVPlacement(
00093 (G4RotationMatrix*)0,
00094 (G4ThreeVector)0,
00095 itsInnerLogVol,
00096 itsName+"_combined",
00097 itsSolidLogVol,
00098 false,
00099 0, BDSGlobalConstants::Instance()->GetCheckOverlaps());
00100
00101 if(BDSGlobalConstants::Instance()->GetSensitiveComponents()){
00102 SetSensitiveVolume(itsSolidLogVol);
00103 }
00104
00105 #ifndef NOUSERLIMITS
00106 itsSolidLogVol->
00107 SetUserLimits(new G4UserLimits(DBL_MAX,DBL_MAX,BDSGlobalConstants::Instance()->GetMaxTime(),
00108 BDSGlobalConstants::Instance()-> GetThresholdCutCharged()));
00109 #endif
00110 itsPhysiComp =
00111 new G4PVPlacement(
00112 (G4RotationMatrix*)0,
00113 (G4ThreeVector)0,
00114 itsSolidLogVol,
00115 itsName+"_solid",
00116 itsMarkerLogicalVolume,
00117 false,
00118 0, BDSGlobalConstants::Instance()->GetCheckOverlaps());
00119 }
00120
00121
00122 BDSSpoiler::~BDSSpoiler()
00123 {
00124 delete itsVisAttributes;
00125
00126
00127
00128
00129
00130
00131 }