/scratch0/jsnuveri/BDSIM/BDSIMgit/bdsim/src/BDSBeamPipeFactoryBase.cc

00001 #include "BDSBeamPipeFactoryBase.hh"
00002 
00003 #include "BDSDebug.hh"
00004 #include "BDSExecOptions.hh"
00005 #include "BDSGlobalConstants.hh"
00006 #include "BDSMaterials.hh"
00007 #include "BDSSDManager.hh"
00008 #include "BDSUtilities.hh"            // for calculateorientation
00009 
00010 #include "globals.hh"                 // geant4 globals / types
00011 #include "G4Colour.hh"
00012 #include "G4LogicalVolume.hh"
00013 #include "G4Material.hh"
00014 #include "G4PVPlacement.hh"
00015 #include "G4ThreeVector.hh"
00016 #include "G4UserLimits.hh"
00017 #include "G4VisAttributes.hh"
00018 
00019 BDSBeamPipeFactoryBase::BDSBeamPipeFactoryBase()
00020 {
00021   lengthSafety              = BDSGlobalConstants::Instance()->GetLengthSafety();
00022   checkOverlaps             = BDSGlobalConstants::Instance()->GetCheckOverlaps();
00023   maxStepFactor             = 0.5; // fraction of length for maximum step size
00024   nSegmentsPerCircle        = 50;
00025   CleanUp();
00026 }
00027 
00028 void BDSBeamPipeFactoryBase::CleanUp()
00029 {
00030   vacuumSolid               = NULL;
00031   beamPipeSolid             = NULL;
00032   containerSolid            = NULL;
00033   containerSubtractionSolid = NULL;
00034   vacuumLV                  = NULL;
00035   beamPipeLV                = NULL;
00036   containerLV               = NULL;
00037 }
00038 
00039 
00040 BDSBeamPipe* BDSBeamPipeFactoryBase::CreateBeamPipeAngledIn(G4String    nameIn,
00041                                                             G4double    lengthIn,
00042                                                             G4double    angleInIn, // the normal angle of the input face
00043                                                             G4double    aper1,
00044                                                             G4double    aper2,
00045                                                             G4double    aper3,
00046                                                             G4double    aper4,
00047                                                             G4Material* vacuumMaterialIn,
00048                                                             G4double    beamPipeThicknessIn,
00049                                                             G4Material* beamPipeMaterialIn)
00050 {
00051   return CreateBeamPipeAngledInOut(nameIn,lengthIn,angleInIn,0,aper1,aper2,aper3,aper4,vacuumMaterialIn,beamPipeThicknessIn,beamPipeMaterialIn);
00052 }
00053 
00054 BDSBeamPipe* BDSBeamPipeFactoryBase::CreateBeamPipeAngledOut(G4String    nameIn,
00055                                                              G4double    lengthIn,
00056                                                              G4double    angleOutIn, // the normal angle of the output face
00057                                                              G4double    aper1,
00058                                                              G4double    aper2,
00059                                                              G4double    aper3,
00060                                                              G4double    aper4,
00061                                                              G4Material* vacuumMaterialIn,
00062                                                              G4double    beamPipeThicknessIn,
00063                                                              G4Material* beamPipeMaterialIn)
00064 {
00065   return CreateBeamPipeAngledInOut(nameIn,lengthIn,0,angleOutIn,aper1,aper2,aper3,aper4,vacuumMaterialIn,beamPipeThicknessIn,beamPipeMaterialIn);
00066 }
00067 
00068 void BDSBeamPipeFactoryBase::TestInputParameters(G4Material*&  vacuumMaterialIn,
00069                                                  G4double&     beamPipeThicknessIn,
00070                                                  G4Material*&  beamPipeMaterialIn)
00071 {
00072   if (!vacuumMaterialIn)
00073     {vacuumMaterialIn = BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->GetVacuumMaterial());}
00074 
00075   if (beamPipeThicknessIn < 1e-10)
00076     {beamPipeThicknessIn = BDSGlobalConstants::Instance()->GetBeamPipeThickness();}
00077 
00078   if (!beamPipeMaterialIn)
00079     {beamPipeMaterialIn = BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->GetBeamPipeMaterialName());}
00080 }
00081   
00082 void BDSBeamPipeFactoryBase::CommonConstruction(G4String    nameIn,
00083                                                 G4Material* vacuumMaterialIn,
00084                                                 G4Material* beamPipeMaterialIn,
00085                                                 G4double    lengthIn)
00086 {
00087 #ifdef BDSDEBUG
00088   G4cout << __METHOD_NAME__ << G4endl;
00089 #endif
00091   BuildLogicalVolumes(nameIn,vacuumMaterialIn,beamPipeMaterialIn);
00093   SetVisAttributes();
00094 #ifndef NOUSERLIMITS
00096   SetUserLimits(lengthIn);
00097 #endif
00099   PlaceComponents(nameIn);
00100 }
00101 
00102 void BDSBeamPipeFactoryBase::BuildLogicalVolumes(G4String    nameIn,
00103                                                  G4Material* vacuumMaterialIn,
00104                                                  G4Material* beamPipeMaterialIn)
00105 {
00106 #ifdef BDSDEBUG
00107   G4cout << __METHOD_NAME__ << G4endl;
00108 #endif
00109   // build the logical volumes
00110   vacuumLV   = new G4LogicalVolume(vacuumSolid,
00111                                    vacuumMaterialIn,
00112                                    nameIn + "_vacuum_lv");
00113   
00114   beamPipeLV = new G4LogicalVolume(beamPipeSolid,
00115                                    beamPipeMaterialIn,
00116                                    nameIn + "_beampipe_lv");
00117   
00118   containerLV = new G4LogicalVolume(containerSolid,
00119                                     vacuumMaterialIn,
00120                                     nameIn + "_container_lv");
00121 }
00122 
00123 void BDSBeamPipeFactoryBase::SetVisAttributes()
00124 {
00125 #ifdef BDSDEBUG
00126   G4cout << __METHOD_NAME__ << G4endl;
00127 #endif
00128   // VISUAL ATTRIBUTES
00129   // set visual attributes
00130   // beampipe
00131   G4VisAttributes* pipeVisAttr = new G4VisAttributes(G4Color(0.4,0.4,0.4));
00132   pipeVisAttr->SetVisibility(true);
00133   pipeVisAttr->SetForceLineSegmentsPerCircle(nSegmentsPerCircle);
00134   beamPipeLV->SetVisAttributes(pipeVisAttr);
00135   // vacuum
00136   vacuumLV->SetVisAttributes(BDSGlobalConstants::Instance()->GetInvisibleVisAttr());
00137   // container
00138   if (BDSExecOptions::Instance()->GetVisDebug())
00139     {containerLV->SetVisAttributes(BDSGlobalConstants::Instance()->GetVisibleDebugVisAttr());}
00140   else
00141     {containerLV->SetVisAttributes(BDSGlobalConstants::Instance()->GetInvisibleVisAttr());}
00142 }
00143 
00144 G4UserLimits* BDSBeamPipeFactoryBase::SetUserLimits(G4double lengthIn)
00145 {
00146 #ifdef BDSDEBUG
00147   G4cout << __METHOD_NAME__ << G4endl;
00148 #endif
00149   // USER LIMITS
00150   // set user limits based on bdsim user specified parameters
00151   G4UserLimits* beamPipeUserLimits = new G4UserLimits("beampipe_cuts");
00152   beamPipeUserLimits->SetMaxAllowedStep( lengthIn * maxStepFactor );
00153   beamPipeUserLimits->SetUserMinEkine(BDSGlobalConstants::Instance()->GetThresholdCutCharged());
00154   beamPipeUserLimits->SetUserMaxTime(BDSGlobalConstants::Instance()->GetMaxTime());
00155   //attach cuts to volumes
00156   vacuumLV->SetUserLimits(beamPipeUserLimits);
00157   beamPipeLV->SetUserLimits(beamPipeUserLimits);
00158   containerLV->SetUserLimits(beamPipeUserLimits);
00159   return beamPipeUserLimits;
00160 }
00161 
00162 void BDSBeamPipeFactoryBase::PlaceComponents(G4String nameIn)
00163 {
00164 #ifdef BDSDEBUG
00165   G4cout << __METHOD_NAME__ << G4endl;
00166 #endif
00167   // PLACEMENT
00168   // place the components inside the container
00169   // note we don't need the pointer for anything - it's registered upon construction with g4
00170   
00171   new G4PVPlacement((G4RotationMatrix*)0,         // no rotation
00172                     (G4ThreeVector)0,             // position
00173                     vacuumLV,                     // lv to be placed
00174                     nameIn + "_vacuum_pv",        // name
00175                     containerLV,                  // mother lv to be place in
00176                     false,                        // no boolean operation
00177                     0,                            // copy number
00178                     checkOverlaps);               // whether to check overlaps
00179   
00180   new G4PVPlacement((G4RotationMatrix*)0,         // no rotation
00181                     (G4ThreeVector)0,             // position
00182                     beamPipeLV,                   // lv to be placed
00183                     nameIn + "_beampipe_pv",      // name
00184                     containerLV,                  // mother lv to be place in
00185                     false,                        // no boolean operation
00186                     0,                            // copy number
00187                     checkOverlaps);               // whether to check overlaps
00188 }
00189 
00190 BDSBeamPipe* BDSBeamPipeFactoryBase::BuildBeamPipeAndRegisterVolumes(std::pair<double,double> extX,
00191                                                                      std::pair<double,double> extY,
00192                                                                      std::pair<double,double> extZ,
00193                                                                      G4double containerRadius)
00194 {  
00195   // build the BDSBeamPipe instance and return it
00196   BDSBeamPipe* aPipe = new BDSBeamPipe(containerSolid,containerLV,extX,extY,extZ,
00197                                        containerSubtractionSolid,
00198                                        vacuumLV,false,containerRadius);
00199 
00200   // REGISTER all lvs
00201   aPipe->RegisterLogicalVolume(vacuumLV); //using geometry component base class method
00202   aPipe->RegisterLogicalVolume(beamPipeLV);
00203 
00204   // register sensitive volumes
00205   aPipe->RegisterSensitiveVolume(beamPipeLV);
00206   
00207   return aPipe;
00208 }
00209 
00211 std::pair<G4ThreeVector,G4ThreeVector> BDSBeamPipeFactoryBase::CalculateFaces(G4double angleIn,
00212                                                                               G4double angleOut)
00213 {
00216   G4int orientationIn  = BDS::CalculateOrientation(angleIn);
00217   G4int orientationOut = BDS::CalculateOrientation(angleOut);
00218   
00219   G4double in_z  = cos(fabs(angleIn)); // calculate components of normal vectors (in the end mag(normal) = 1)
00220   G4double in_x  = sin(fabs(angleIn)); // note full angle here as it's the exit angle
00221   G4double out_z = cos(fabs(angleOut));
00222   G4double out_x = sin(fabs(angleOut));
00223   G4ThreeVector inputface  = G4ThreeVector(orientationIn*in_x, 0.0, -1.0*in_z); //-1 as pointing down in z for normal
00224   G4ThreeVector outputface = G4ThreeVector(orientationOut*out_x, 0.0, out_z);   // no output face angle
00225   return std::make_pair(inputface,outputface);
00226 }

Generated on 28 Jun 2015 for BDSIM by  doxygen 1.4.7