include/BDSAcceleratorComponent.hh

00001 //  
00002 //   BDSIM, (C) 2001-2006 
00003 //   
00004 //   version 0.3
00005 //  
00006 //
00007 //
00008 //
00009 //
00010 //   Generic accelerator component class
00011 //
00012 //
00013 //   History
00014 //
00015 //     24 Nov 2006 by Agapov,  v.0.3
00016 //     x  x   2002 by Blair
00017 //
00018 //
00019 
00020 
00021 
00022 
00023 #ifndef __BDSACCELERATORCOMPONENT_H
00024 #define __BDSACCELERATORCOMPONENT_H
00025 
00026 #include "BDSGlobalConstants.hh" 
00027 
00028 #include <cstring>
00029 #include <list>
00030 #include <vector> 
00031 #include "G4LogicalVolume.hh"
00032 #include "G4VisAttributes.hh"
00033 #include "G4LogicalVolume.hh"
00034 #include "globals.hh"
00035 //#include "BDSBeamPipe.hh"
00036 #include "BDSEnergyCounterSD.hh"
00037 
00038 #include "G4MagneticField.hh"
00039 #include "G4Mag_EqRhs.hh"
00040 #include "G4MagIntegratorStepper.hh"
00041 #include "G4FieldManager.hh"
00042 #include "G4UserLimits.hh"
00043 #include "G4PVPlacement.hh"
00044 #include "G4AssemblyVolume.hh"
00045 #include "G4CSGSolid.hh"
00046 #include "G4Tubs.hh"
00047 
00048 
00049 class BDSAcceleratorComponent 
00050 {
00051 public:
00052   //destructor
00053   virtual ~BDSAcceleratorComponent ();
00054 
00055   //name
00056   const G4String GetName () const;
00057   void SetName(G4String aName);
00058 
00059   //type 
00060   const G4String GetType () const;
00061   void SetType(G4String aType);
00062 
00063   G4int GetPrecisionRegion() const; //0 = no precision region, 1 = precision region 1, 2 = precision region 2.
00064   void SetPrecisionRegion(G4int aPrecisionRegion);
00065   
00066 
00067   //
00068   //    Geometry features    
00069   //
00070 
00071   void BuildTunnel();
00072   void BuildBLMs();
00073   void BuildGate();
00074 
00075   // angle - for bends etc.
00076   G4double GetAngle ();
00077 
00078   // geometry length of the component.
00079   void SetLength(G4double aLength); 
00080   virtual G4double GetLength ();
00081   virtual G4double GetZLength ();
00082   virtual G4double GetXLength ();
00083   virtual G4double GetYLength ();
00084   virtual G4double GetArcLength ();
00085 
00086   G4double GetPhiAngleIn (); //polar angle in
00087   G4double GetPhiAngleOut (); //polar angle out
00088 
00089   G4double GetPhi (); //polar angle with respect to original frame
00090   void SetPhi (G4double val);
00091 
00092   G4double GetTheta (); //azimuthal angle with respect to original frame
00093   void SetTheta(G4double val);
00094 
00095   G4double GetPsi (); //azimuthal angle with respect to original frame
00096   void SetPsi(G4double val);
00097 
00098   G4double GetXOffset();  // frame offset 
00099   G4double GetYOffset();
00100   G4double GetZOffset();
00101 
00102   G4double GetTunnelRadius();
00103   G4double GetTunnelOffsetX();
00104   
00105   G4double GetAperX();
00106   G4double GetAperY();
00107 
00108   G4double GetK1();
00109   G4double GetK2();
00110   G4double GetK3();
00111 
00112   //Set is only for Outline readout purposes - doesn't change magnet strengths
00113   void SetK1(G4double K1);
00114   void SetK2(G4double K2);
00115   void SetK3(G4double K3);
00116 
00117   G4RotationMatrix* GetRotation();
00118   G4ThreeVector GetPosition();
00119   void SetPosition(G4ThreeVector);
00120   
00121   G4double GetTilt();  // component tilt 
00122   
00123   
00124   
00125 
00126 
00127   G4LogicalVolume* GetMarkerLogicalVolume() const;
00128 
00129   G4LogicalVolume* GetTunnelLogicalVolume() const;
00130   G4String GetTunnelCavityMaterial() const;
00131   BDSEnergyCounterSD* GetBDSEnergyCounter() const;
00132   
00133   void SetBDSEnergyCounter( BDSEnergyCounterSD* anBDSEnergyCounter);
00134   
00135   G4int GetCopyNumber() const;
00136   
00137   G4double GetSPos() const;
00138   
00139   void SetSPos(G4double spos);
00140   
00141   void SetCopyNumber(G4int nCopy);
00142   
00143   void SetSensitiveVolume(G4LogicalVolume* aLogVol);
00144   
00145   G4LogicalVolume* GetSensitiveVolume();
00146 
00147   void SetMultipleSensitiveVolumes(G4LogicalVolume* aLogVol);
00148 
00149   std::vector<G4LogicalVolume*> GetMultipleSensitiveVolumes();
00150 
00151   void SetGFlashVolumes(G4LogicalVolume* aLogVol);
00152 
00153   std::vector<G4LogicalVolume*> GetGFlashVolumes();
00154 
00155   void SetMultiplePhysicalVolumes(G4VPhysicalVolume* aPhysVol);
00156 
00157   std::vector<G4VPhysicalVolume*> GetMultiplePhysicalVolumes();
00158 
00159   void SetInnerMostLogicalVolume(G4LogicalVolume* aLogVol);
00160   
00161   G4LogicalVolume* GetInnerMostLogicalVolume() const;
00162   
00163   G4UserLimits* GetInnerBPUserLimits();
00164   G4UserLimits* GetUserLimits();
00165 
00166   // G4double GetZLower();
00167   // G4double GetZUpper();
00168   // void SetZLower(G4double aZLower);
00169   // void SetZUpper(G4double aZUpper);
00170   // void AddSynchEnergyLoss(G4double SynchEnergyLoss);
00171   // G4double GetSynchEnergyLoss();
00172   
00173   BDSAcceleratorComponent();
00174   void BuildOuterFieldManager();
00175 
00176   // in case a mapped field is provided creates a field mesh in global coordinates
00177   virtual void PrepareField(G4VPhysicalVolume *referenceVolume); 
00178   virtual void SynchRescale(G4double factor); 
00179 
00180   // in case a component requires specific alignment (e.g. SQL/BDSElement)
00181   virtual void AlignComponent(G4ThreeVector& TargetPos, 
00182                               G4RotationMatrix *TargetRot,
00183                               G4RotationMatrix& globalRotation,
00184                               G4ThreeVector& rtot,
00185                               G4ThreeVector& rlast,
00186                               G4ThreeVector& localX,
00187                               G4ThreeVector& localY,
00188                               G4ThreeVector& localZ); 
00189 
00190   
00191   // get parameter value from the specification string
00192 
00193   G4double getParameterValue(G4String spec, G4String name) const;
00194   G4String getParameterValueString(G4String spec, G4String name) const;
00195 
00196   // constructor
00197   BDSAcceleratorComponent (
00198                           G4String& aName, 
00199                           G4double aLength,
00200                           G4double aBpRadius,
00201                           G4double aXAper,
00202                           G4double aYAper,
00203                           G4VisAttributes* aVisAtt,
00204                           std::list<G4double> blmLocZ, std::list<G4double> blmLocTheta,
00205                           G4String aTunnelMaterial = "",
00206                           G4String aMaterial = "",
00207                           G4double phi=0.,  // polar angle (used in hor. bends)
00208                           //G4double theta=0.,
00209                           G4double XOffset=0.,
00210                           G4double YOffset=0.,
00211                           G4double ZOffset=0.,
00212                           G4double tunnelRadius=0.,
00213                           G4double tunnelOffsetX=BDSGlobalConstants::Instance()->GetTunnelOffsetX(),
00214                           G4String aTunnelCavityMaterial = "Air",
00215                           G4int aPrecisionRegion=0
00216 );
00217 
00218   BDSAcceleratorComponent (
00219                           G4String& aName, 
00220                           G4double aLength,
00221                           G4double aBpRadius,
00222                           G4double aXAper,
00223                           G4double aYAper,
00224                           G4VisAttributes* aVisAtt,
00225                           G4String aTunnelMaterial = "",
00226                           G4String aMaterial = "",
00227                           G4double phi=0.,  // polar angle (used in hor. bends)
00228                           G4double XOffset=0.,
00229                           G4double YOffset=0.,
00230                           G4double ZOffset=0.,
00231                           G4double tunnelRadius=0.,
00232                           G4double tunnelOffsetX=BDSGlobalConstants::Instance()->GetTunnelOffsetX(),
00233                           G4String aTunnelCavityMaterial = "Air",
00234                           G4int aPrecisionRegion=0);
00235 
00236 
00237 
00238   G4VisAttributes* GetVisAttributes()const;
00239   G4LogicalVolume* itsOuterLogicalVolume;
00240   G4LogicalVolume* itsMarkerLogicalVolume;
00241   G4LogicalVolume* itsTunnelLogicalVolume;
00242   G4LogicalVolume* itsTunnelFloorLogicalVolume;
00243 
00244 
00245 protected:
00246 
00247   //Calculate dimensions used for the marker volume etc.
00248   void CalculateLengths();
00249 
00250 
00251   //Values related to BLM placement and geometry
00252   G4double itsBlmLocationR;
00253   //  G4double itsBlmRadius;
00254 
00255   G4String itsName;
00256   G4double itsLength;
00257   G4double itsXLength;
00258   G4double itsYLength;
00259   G4double itsOuterR;
00260   G4double itsBpRadius;
00261   G4double itsXAper;
00262   G4double itsYAper;
00263   G4double itsAngle;
00264   G4String itsMaterial;
00265   G4VisAttributes* itsVisAttributes;
00266   std::list<G4double> itsBlmLocZ;
00267   std::list<G4double> itsBlmLocTheta;
00268   G4String itsTunnelMaterial;
00269   //Tunnel geom
00270   G4double itsXOffset;
00271   G4double itsYOffset;
00272   G4double itsZOffset;
00273   G4double itsTunnelRadius;
00274   G4double itsTunnelOffsetX;  
00275 
00276   G4String itsType;
00277 
00278   G4double itsTilt;
00279 
00280   G4double itsPhiAngleIn;
00281   G4double itsPhiAngleOut;
00282   
00283   G4double itsMagScaleFactor;
00284   G4double itsPhi;
00285   G4double itsTheta;
00286   G4double itsPsi;
00287   G4double itsK1, itsK2, itsK3;
00288   G4RotationMatrix* itsRotation;
00289   G4ThreeVector itsPosition;
00290   //  BDSBeamPipe* itsBeamPipe;
00291   G4MagIntegratorStepper*  itsOuterStepper;
00293   G4UserLimits* itsUserLimits;
00295   G4UserLimits* itsOuterUserLimits;
00296   G4UserLimits* itsMarkerUserLimits;
00297   G4UserLimits* itsInnerBeampipeUserLimits;
00298   G4LogicalVolume* itsInnerMostLogicalVolume;
00299 
00300   G4String itsTunnelCavityMaterial;
00301   G4int itsPrecisionRegion;
00302 
00304   G4VSolid* itsMarkerSolidVolume;
00305 
00306 
00307 
00309   G4VSolid* itsTunnelSolid;
00310   G4VSolid* itsSoilSolid;
00311   G4VSolid* itsInnerTunnelSolid;
00312   G4VSolid *itsTunnelCavity;
00313   G4VSolid *itsLargerTunnelCavity;
00314   G4VSolid *itsTunnelFloor;
00315   G4VSolid* itsLargerInnerTunnelSolid; 
00316   G4VSolid *itsTunnelMinusCavity;
00317   G4CSGSolid* itsTunnelSizedBlock;
00318 
00320   G4LogicalVolume* itsBLMLogicalVolume;
00321   G4LogicalVolume* itsBlmCaseLogicalVolume;
00323   std::vector<G4VPhysicalVolume*> itsBLMPhysiComp;
00325   G4LogicalVolume* itsSoilTunnelLogicalVolume;
00326   G4LogicalVolume* itsTunnelCavityLogicalVolume;
00327   G4LogicalVolume*  itsTunnelMinusCavityLogicalVolume;
00329   G4VPhysicalVolume* itsTunnelPhysiInner;
00330   G4VPhysicalVolume* itsTunnelPhysiComp;
00331   G4VPhysicalVolume* itsTunnelFloorPhysiComp;
00332   G4VPhysicalVolume* itsTunnelPhysiCompSoil;
00334   G4UserLimits* itsTunnelUserLimits;
00335   G4UserLimits* itsSoilTunnelUserLimits;
00336   G4UserLimits* itsInnerTunnelUserLimits;
00337 
00338 
00339 
00340 private:
00342   void ConstructorInit();
00343 
00344   G4RotationMatrix* nullRotationMatrix;
00345   G4RotationMatrix* tunnelRot;
00346   G4VisAttributes* VisAtt;
00347   G4VisAttributes* VisAtt1;
00348   G4VisAttributes* VisAtt2;
00349   G4VisAttributes* VisAtt3;
00350   G4VisAttributes* VisAtt4;
00351   G4VisAttributes* VisAtt5;
00352   G4Tubs* itsBLMSolid;
00353   G4Tubs* itsBlmOuterSolid;
00354   G4double itsSPos;
00355   G4int itsCopyNumber;
00356   BDSEnergyCounterSD* itsBDSEnergyCounter;
00357   //  G4int itsCollectionID;
00358   G4LogicalVolume* itsSensitiveVolume;
00359   std::vector<G4LogicalVolume*> itsMultipleSensitiveVolumes;
00360   std::vector<G4LogicalVolume*> itsGFlashVolumes;
00361   //A vector containing the physical volumes in the accelerator component- to be used for geometric importance sampling etc.
00362   std::vector<G4VPhysicalVolume*> itsMultiplePhysicalVolumes;
00363   //  G4double itsZLower;
00364   //  G4double itsZUpper;
00365   //  G4double itsSynchEnergyLoss;
00366 
00367 };
00368 
00369 // Class BDSAcceleratorComponent 
00370 
00371 inline G4double BDSAcceleratorComponent::GetLength ()
00372 {return itsLength;}
00373 
00374 inline G4double BDSAcceleratorComponent::GetXLength ()
00375 {return itsXLength;}
00376 
00377 inline G4double BDSAcceleratorComponent::GetYLength ()
00378 {return itsYLength;}
00379 
00380 inline G4double BDSAcceleratorComponent::GetArcLength ()
00381 {return itsLength;}
00382 
00383 inline G4double BDSAcceleratorComponent::GetZLength ()
00384 {return itsLength;}
00385 
00386 inline G4double BDSAcceleratorComponent::GetAngle ()
00387 {return itsAngle;}
00388 
00389 inline G4double BDSAcceleratorComponent::GetPhiAngleIn ()
00390 {return itsPhiAngleIn;}
00391 
00392 inline G4double BDSAcceleratorComponent::GetPhiAngleOut ()
00393 {return itsPhiAngleOut;}
00394 
00395 inline G4double BDSAcceleratorComponent::GetPhi ()
00396 {return itsPhi;}
00397 
00398 inline void BDSAcceleratorComponent::SetPhi (G4double val)
00399 {itsPhi = val;}
00400 
00401 inline G4double BDSAcceleratorComponent::GetTheta ()
00402 {return itsTheta;}
00403 
00404 inline void BDSAcceleratorComponent::SetTheta (G4double val)
00405 {itsTheta = val;}
00406 
00407 inline G4double BDSAcceleratorComponent::GetPsi ()
00408 {return itsPsi;}
00409 
00410 inline void BDSAcceleratorComponent::SetPsi (G4double val)
00411 {itsPsi = val;}
00412 
00413 inline G4double BDSAcceleratorComponent::GetAperX()
00414 {
00415   if(itsXAper==0) // i.e. it has not been set
00416     return itsBpRadius;
00417   else return itsXAper;
00418 }
00419 
00420 inline G4double BDSAcceleratorComponent::GetAperY()
00421 {
00422   if(itsYAper==0) // i.e. it has not been set
00423     return itsBpRadius;
00424   else return itsYAper;
00425 }
00426 
00427 inline G4double BDSAcceleratorComponent::GetK1()
00428 { return itsK1; }
00429 
00430 inline G4double BDSAcceleratorComponent::GetK2()
00431 { return itsK2; }
00432 
00433 inline G4double BDSAcceleratorComponent::GetK3()
00434 { return itsK3; }
00435 
00436 inline void BDSAcceleratorComponent::SetK1(G4double K1)
00437 { itsK1 = K1; }
00438 
00439 inline void BDSAcceleratorComponent::SetK2(G4double K2)
00440 { itsK2 = K2; }
00441 
00442 inline void BDSAcceleratorComponent::SetK3(G4double K3)
00443 { itsK3 = K3; }
00444 
00445 inline G4RotationMatrix* BDSAcceleratorComponent::GetRotation()
00446 { return itsRotation;}
00447 
00448 inline G4ThreeVector BDSAcceleratorComponent::GetPosition()
00449 { return itsPosition;}
00450 
00451 inline void BDSAcceleratorComponent::SetPosition(G4ThreeVector pos)
00452 { itsPosition = pos;}
00453 
00454 inline const G4String BDSAcceleratorComponent::GetName () const
00455 {return itsName;}
00456 
00457 inline void BDSAcceleratorComponent::SetName (G4String aName)
00458 {itsName=aName;}
00459 
00460 inline const G4String BDSAcceleratorComponent::GetType () const
00461 {return itsType;}
00462 
00463 inline void BDSAcceleratorComponent::SetType (G4String aType)
00464 {itsType=aType;}
00465 
00466 inline G4int BDSAcceleratorComponent::GetPrecisionRegion () const
00467 {return itsPrecisionRegion;}
00468 
00469 inline void BDSAcceleratorComponent::SetPrecisionRegion (G4int aPrecisionRegion)
00470 {itsPrecisionRegion=aPrecisionRegion;}
00471 
00472 inline G4LogicalVolume* BDSAcceleratorComponent::GetMarkerLogicalVolume() const
00473 {return itsMarkerLogicalVolume;}
00474 
00475 inline G4LogicalVolume* BDSAcceleratorComponent::GetInnerMostLogicalVolume() const
00476 {return itsInnerMostLogicalVolume;}
00477 
00478 inline void BDSAcceleratorComponent::
00479 SetInnerMostLogicalVolume(G4LogicalVolume* aLogVol)
00480 {itsInnerMostLogicalVolume = aLogVol;}
00481 
00482 inline G4VisAttributes* BDSAcceleratorComponent::GetVisAttributes() const
00483 {return itsVisAttributes;}
00484 
00485 inline BDSEnergyCounterSD* BDSAcceleratorComponent::GetBDSEnergyCounter() const
00486 {return itsBDSEnergyCounter;}
00487 
00488 inline G4int BDSAcceleratorComponent::GetCopyNumber() const
00489 {return itsCopyNumber;}
00490 
00491 inline G4double BDSAcceleratorComponent::GetSPos() const
00492 {return itsSPos;}
00493 
00494 inline void BDSAcceleratorComponent::SetCopyNumber(G4int nCopy)
00495 {itsCopyNumber=nCopy;}
00496 
00497 inline void BDSAcceleratorComponent::SetSPos(G4double spos)
00498 {itsSPos=spos;}
00499 
00500 inline void 
00501 BDSAcceleratorComponent::SetBDSEnergyCounter(BDSEnergyCounterSD* anBDSEnergyCounter)
00502 {itsBDSEnergyCounter=anBDSEnergyCounter;}
00503 
00504 inline  
00505 void BDSAcceleratorComponent::SetSensitiveVolume(G4LogicalVolume* aLogVol)
00506 {itsSensitiveVolume=aLogVol;}
00507 
00508 inline  G4LogicalVolume* BDSAcceleratorComponent::GetSensitiveVolume()
00509 {return itsSensitiveVolume;}
00510 
00511 inline void BDSAcceleratorComponent::SetMultipleSensitiveVolumes(G4LogicalVolume* aLogVol)
00512 { itsMultipleSensitiveVolumes.push_back(aLogVol);}
00513 
00514 inline void BDSAcceleratorComponent::SetGFlashVolumes(G4LogicalVolume* aLogVol)
00515 { itsGFlashVolumes.push_back(aLogVol);}
00516 
00517 inline  std::vector<G4LogicalVolume*> BDSAcceleratorComponent::GetMultipleSensitiveVolumes()
00518 {return itsMultipleSensitiveVolumes;}
00519 
00520 inline  std::vector<G4LogicalVolume*> BDSAcceleratorComponent::GetGFlashVolumes()
00521 {return itsGFlashVolumes;}
00522 
00523 inline void BDSAcceleratorComponent::SetMultiplePhysicalVolumes(G4VPhysicalVolume* aPhysVol)
00524 { itsMultiplePhysicalVolumes.push_back(aPhysVol);}
00525 
00526 inline  std::vector<G4VPhysicalVolume*> BDSAcceleratorComponent::GetMultiplePhysicalVolumes()
00527 {return itsMultiplePhysicalVolumes;}
00528 
00529 // inline  G4double BDSAcceleratorComponent::GetZLower()
00530 // {return itsZLower;}
00531 
00532 // inline  G4double BDSAcceleratorComponent::GetZUpper()
00533 // {return itsZUpper;}
00534 
00535 // inline  void BDSAcceleratorComponent::SetZLower(G4double aZLower)
00536 // {itsZLower=aZLower;}
00537 
00538 // inline  void BDSAcceleratorComponent::SetZUpper(G4double aZUpper)
00539 // {itsZUpper=aZUpper;}
00540 
00541 // inline void 
00542 // BDSAcceleratorComponent::AddSynchEnergyLoss(G4double SynchEnergyLoss)
00543 // {itsSynchEnergyLoss+=SynchEnergyLoss;}
00544 
00545 // inline  G4double BDSAcceleratorComponent::GetSynchEnergyLoss()
00546 // {return itsSynchEnergyLoss;}
00547 
00548 inline  G4UserLimits* BDSAcceleratorComponent::GetUserLimits(){
00549   return itsUserLimits;
00550 }
00551 
00552 inline  G4UserLimits* BDSAcceleratorComponent::GetInnerBPUserLimits()
00553   {return itsInnerBeampipeUserLimits;}
00554 
00555 inline  G4double BDSAcceleratorComponent::GetXOffset()
00556 {return itsXOffset;}
00557 
00558 inline G4double BDSAcceleratorComponent::GetYOffset() 
00559 {return itsYOffset;}
00560 
00561 inline G4double BDSAcceleratorComponent::GetZOffset()
00562 {return itsZOffset;}
00563 
00564 inline G4double BDSAcceleratorComponent::GetTunnelRadius()
00565 {return itsTunnelRadius;}
00566 
00567 inline G4double BDSAcceleratorComponent::GetTunnelOffsetX()
00568 {return itsTunnelOffsetX;}
00569 
00570 inline G4double BDSAcceleratorComponent::GetTilt()
00571 {return itsTilt;}
00572 
00573 
00574 inline  G4double BDSAcceleratorComponent::getParameterValue(G4String spec, G4String name) const
00575 {
00576   G4double value = 0;
00577 
00578   std::string delimiters = "&";
00579   std::string param = name + "=";
00580 
00581   int pos = spec.find(param);
00582   if( pos >= 0 )
00583     {
00584       
00585       int pos2 = spec.find("&",pos);
00586       int pos3 = spec.length();
00587       int tend = pos2 < 0 ? pos3 : pos2; 
00588       int llen = tend - pos - param.length();
00589       
00590       std::string val = spec.substr(pos + param.length(), llen);
00591       
00592       value = atof(val.c_str());
00593 
00594   }
00595 
00596   return value;
00597 
00598 }
00599 
00600 inline  G4String BDSAcceleratorComponent::getParameterValueString(G4String spec, G4String name) const
00601 {
00602   G4String value = "";
00603 
00604   std::string delimiters = "&";
00605   std::string param = name + "=";
00606 
00607   int pos = spec.find(param);
00608   if( pos >= 0 )
00609     {
00610       
00611       int pos2 = spec.find("&",pos);
00612       int pos3 = spec.length();
00613       int tend = pos2 < 0 ? pos3 : pos2; 
00614       int llen = tend - pos - param.length();
00615       
00616       value = spec.substr(pos + param.length(), llen);
00617   }
00618 
00619   return value;
00620 
00621 }
00622 
00623 #endif

Generated on 27 Aug 2013 for BDSIM by  doxygen 1.4.7