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

00001 #include "BDSGlobalConstants.hh" 
00002 #include "BDSCCDCamera.hh"
00003 #include "BDSCCDPixel.hh"
00004 #include "G4VisAttributes.hh"
00005 #include "G4LogicalVolume.hh"
00006 #include "G4VPhysicalVolume.hh"
00007 #include "G4PVPlacement.hh"               
00008 #include "G4Box.hh"
00009 #include "G4Tubs.hh"
00010 #include "BDSDebug.hh"
00011 
00012 #include <map>
00013 #include "G4TwoVector.hh"
00014 #include "BDSMaterials.hh"
00015 
00016 BDSCCDCamera::BDSCCDCamera ()
00017 {
00018   _name="ccdCamera";
00019   defaultDimensions();
00020   build();
00021 }
00022 
00023 void BDSCCDCamera::defaultDimensions(){
00024   _cavityLength = 40*CLHEP::mm;
00025   _size.setX(110.7*CLHEP::mm);
00026   _size.setY(137.2*CLHEP::mm);
00027   _size.setZ(231.0*CLHEP::mm+_cavityLength);
00028 }
00029 
00030 void BDSCCDCamera::build(){
00031   buildMotherVolume();
00032   buildObjectLens();
00033   buildImageLens();
00034   buildCCDChip();
00035   buildCavity();
00036   placeComponents();
00037 }
00038 
00039 void BDSCCDCamera::buildMotherVolume(){
00040   G4cout << __METHOD_NAME__ << G4endl;
00041   _solid=new G4Box( _name+"_solid",
00042                     _size.x()/2.0,
00043                     _size.y()/2.0,
00044                     _size.z()/2.0);
00045   
00046   _log=new G4LogicalVolume
00047     (_solid, 
00048      BDSMaterials::Instance()->GetMaterial("G4_POLYSTYRENE"),
00049      _name+"_log");
00050   _log->SetVisAttributes(new G4VisAttributes(G4Color(0,1.0,0)));
00051   G4cout << __METHOD_END__ << G4endl;
00052 }
00053 
00054 void BDSCCDCamera::buildCavity(){
00055   //  G4double xSize=std::max(_objectLens->diameter(),_imageLens->diameter());
00056   //  xSize=std::max(xSize,_ccdChip->size().x());
00057   //  G4double ySize=std::max(_objectLens->diameter(),_imageLens->diameter());
00058   //  ySize=std::max(ySize,_ccdChip->size().y());
00059  
00060   G4Tubs* cavityTubs = new G4Tubs("CCDCameraCavity_solid",0,_objectLens->diameter()/2.0,_cavityLength/2.0,0,CLHEP::twopi*CLHEP::rad);
00061   _cavityLog = new G4LogicalVolume(                                
00062                                    cavityTubs,     
00063                                    BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->GetVacuumMaterial()),
00064                                    _name+"_cavity_log");
00065   _cavityLog->SetVisAttributes(new G4VisAttributes(G4Color(1.0,0.0,0)));
00066 }
00067 
00068 void BDSCCDCamera::placeCavity(){
00069   G4ThreeVector placementVec;
00070   placementVec.setX(0);
00071   placementVec.setY(0);
00072   placementVec.setZ(-_size.z()/2.0+_cavityLength/2.0);
00073   new G4PVPlacement(0,
00074                     placementVec,
00075                     _cavityLog,
00076                     _name+"_cavity_phys",
00077                     _log,
00078                     false,
00079                     0,
00080                     false);
00081 }
00082 
00083 void BDSCCDCamera::buildObjectLens(){
00084   G4cout << __METHOD_NAME__ << G4endl;
00085   G4double factor =1.0;
00086   _objectLens = new BDSLens(_name+"ObjectLens",factor*25.4*CLHEP::mm, 1029.8*CLHEP::mm,factor*2.2*CLHEP::mm); //Focal length 1m lens (Thorlabs LB1409-A)
00087   G4cout << __METHOD_END__ << G4endl;
00088 }
00089 
00090 void BDSCCDCamera::buildImageLens(){
00091   G4cout << __METHOD_NAME__ << G4endl;
00092   //  _imageLens = new BDSLens(_name+"ImageLens",12.7*CLHEP::mm, 14.6*CLHEP::mm, 4.7*CLHEP::mm); //Focal length 15CLHEP::mm lens (Thorlabs LB1092-A) (magnification factor = 66.4)
00093   _imageLens = new BDSLens(_name+"ImageLens",25.4*CLHEP::mm, 25.5*CLHEP::mm, 9.0*CLHEP::mm); //Focal length 25.4CLHEP::mm lens (Thorlabs LB1761-A)(back focal length 22.2CLHEP::mm)
00094   G4cout << __METHOD_END__ << G4endl;
00095 }
00096 void BDSCCDCamera::buildCCDChip(){
00097   G4cout << __METHOD_NAME__ << G4endl;
00098   G4ThreeVector pixelSize;
00099   G4TwoVector nPixels;
00100 
00101   pixelSize.setX(13.5e-3*CLHEP::mm);
00102   pixelSize.setY(13.5e-3*CLHEP::mm*512);
00103   pixelSize.setZ(13.5e-3*CLHEP::mm); //Assume that the pixels are cubes for the moment.
00104   nPixels.setX(2048);
00105   nPixels.setY(1);
00106 
00107   G4cout << __METHOD_NAME__ << " constructing BDSCCDChip..." << G4endl;
00108   _ccdChip = new BDSCCDChip((G4String)(_name+"_CCDChip"), pixelSize, nPixels);
00109   G4cout << __METHOD_END__ << G4endl;
00110 }
00111 
00112 void BDSCCDCamera::placeComponents(){
00113   placeCavity();
00114   placeObjectLens();
00115   placeImageLens();
00116   placeCCDChip();
00117 }
00118 
00119 void BDSCCDCamera::placeObjectLens(){
00120   G4ThreeVector placementVec;
00121   placementVec.setX(0);
00122   placementVec.setY(0);
00123   placementVec.setZ(-_cavityLength/2.0+_objectLens->centreThickness()/2.0);
00124   new G4PVPlacement(0,
00125                     placementVec,
00126                     _objectLens->log(),
00127                     _objectLens->name()+"_phys",
00128                     _cavityLog,
00129                     false,
00130                     0,
00131                     false);
00132 }
00133 
00134 void BDSCCDCamera::placeImageLens(){
00135   G4cout << __METHOD_NAME__ << G4endl;
00136   G4ThreeVector placementVec;
00137   placementVec.setX(0);
00138   placementVec.setY(0);
00139   placementVec.setZ(_cavityLength/2.0-_imageLens->centreThickness()/2.0-22.2*CLHEP::mm);
00140   new G4PVPlacement(0,
00141                     placementVec,
00142                     _imageLens->log(),
00143                     _imageLens->name()+"_phys",
00144                     _cavityLog,
00145                     false,
00146                     0,
00147                     false);
00148   G4cout << __METHOD_END__ << G4endl;
00149 }
00150 void BDSCCDCamera::placeCCDChip(){
00151   G4cout << __METHOD_NAME__ << G4endl;
00152   G4ThreeVector placementVec;
00153   placementVec.setX(0);
00154   placementVec.setY(0);
00155   placementVec.setZ(_cavityLength/2.0-_ccdChip->size().z()/2.0);
00156   new G4PVPlacement(0,
00157                     placementVec,
00158                     _ccdChip->log(),
00159                     _ccdChip->name()+"_phys",
00160                     _cavityLog,
00161                     false,
00162                     0,
00163                     false);
00164   G4cout << __METHOD_END__ << G4endl;
00165 }
00166 
00167 BDSCCDCamera::~BDSCCDCamera()
00168 {
00169   delete _objectLens;
00170   delete _imageLens;
00171   delete _ccdChip;
00172 
00173 }
00174 

Generated on 28 Jun 2015 for BDSIM by  doxygen 1.4.7