00001 #ifndef __BDSACCELERATORCOMPONENT_H
00002 #define __BDSACCELERATORCOMPONENT_H
00003
00004 #include "G4LogicalVolume.hh"
00005 #include "globals.hh"
00006
00007 #include "BDSBeamPipeInfo.hh"
00008 #include "BDSGeometryComponent.hh"
00009 #include "BDSGlobalConstants.hh"
00010 #include "BDSTiltOffset.hh"
00011
00012 #include <vector>
00013
00042 class BDSAcceleratorComponent: public BDSGeometryComponent
00043 {
00044 public:
00054 BDSAcceleratorComponent(G4String name,
00055 G4double arcLength,
00056 G4double angle,
00057 G4String type,
00058 BDSTiltOffset tiltOffset = BDSTiltOffset(),
00059 G4int precisionRegion = 0,
00060 BDSBeamPipeInfo* beamPipeInfo = NULL);
00061
00062 virtual ~BDSAcceleratorComponent();
00063
00065 G4String GetName() const;
00066
00068 G4String GetType() const;
00069
00071 G4int GetPrecisionRegion() const;
00072
00075 G4double GetAngle() const;
00076
00078 inline BDSTiltOffset GetTiltOffset() const;
00079
00081 virtual G4double GetArcLength() const;
00082 virtual G4double GetChordLength() const;
00083
00085 inline G4LogicalVolume* GetReadOutLogicalVolume() const;
00086
00087
00088 virtual void PrepareField(G4VPhysicalVolume *referenceVolume);
00089
00091 G4double GetParameterValue (G4String spec, G4String name) const;
00092 G4String GetParameterValueString(G4String spec, G4String name) const;
00094
00096 friend class BDSComponentFactory;
00097 friend class BDSLine;
00098 friend class BDSDetectorConstruction;
00099
00101 G4double GetSPos() const;
00102 void SetSPos(G4double spos);
00103 void SetGFlashVolumes(G4LogicalVolume* aLogVol);
00104 std::vector<G4LogicalVolume*> GetGFlashVolumes() const;
00105 void SetMultiplePhysicalVolumes(G4VPhysicalVolume* aPhysVol);
00106 std::vector<G4VPhysicalVolume*> GetMultiplePhysicalVolumes() const;
00108
00109 protected:
00112
00113 virtual void Initialise();
00114
00119 virtual void Build();
00120
00123 virtual void BuildContainerLogicalVolume() = 0;
00124
00126 const G4String name;
00127 const G4double arcLength;
00128 const G4String type;
00130
00132 G4double chordLength;
00133 G4double angle;
00134 BDSTiltOffset tiltOffset;
00135 G4int precisionRegion;
00136 BDSBeamPipeInfo* beamPipeInfo;
00138
00140 static G4double lengthSafety;
00141 static G4Material* emptyMaterial;
00142
00143 private:
00147 BDSAcceleratorComponent();
00148
00150 BDSAcceleratorComponent& operator=(const BDSAcceleratorComponent&);
00151 BDSAcceleratorComponent(BDSAcceleratorComponent&);
00152
00153 G4LogicalVolume* readOutLV;
00154
00156 G4double itsSPos;
00157
00158 std::vector<G4LogicalVolume*> itsGFlashVolumes;
00159
00160
00161
00162 std::vector<G4VPhysicalVolume*> itsMultiplePhysicalVolumes;
00163 };
00164
00165 inline G4String BDSAcceleratorComponent::GetName() const
00166 {return name;}
00167
00168 inline G4double BDSAcceleratorComponent::GetChordLength() const
00169 {return chordLength;}
00170
00171 inline G4double BDSAcceleratorComponent::GetArcLength() const
00172 {return arcLength;}
00173
00174 inline G4double BDSAcceleratorComponent::GetAngle() const
00175 {return angle;}
00176
00177 inline G4String BDSAcceleratorComponent::GetType() const
00178 {return type;}
00179
00180 inline G4int BDSAcceleratorComponent::GetPrecisionRegion() const
00181 {return precisionRegion;}
00182
00183 inline BDSTiltOffset BDSAcceleratorComponent::GetTiltOffset() const
00184 {return tiltOffset;}
00185
00186 inline G4double BDSAcceleratorComponent::GetSPos() const
00187 {return itsSPos;}
00188
00189 inline void BDSAcceleratorComponent::SetSPos(G4double spos)
00190 {itsSPos=spos;}
00191
00192 inline void BDSAcceleratorComponent::SetGFlashVolumes(G4LogicalVolume* aLogVol)
00193 {itsGFlashVolumes.push_back(aLogVol);}
00194
00195 inline std::vector<G4LogicalVolume*> BDSAcceleratorComponent::GetGFlashVolumes() const
00196 {return itsGFlashVolumes;}
00197
00198 inline void BDSAcceleratorComponent::SetMultiplePhysicalVolumes(G4VPhysicalVolume* aPhysVol)
00199 {itsMultiplePhysicalVolumes.push_back(aPhysVol);}
00200
00201 inline std::vector<G4VPhysicalVolume*> BDSAcceleratorComponent::GetMultiplePhysicalVolumes() const
00202 {return itsMultiplePhysicalVolumes;}
00203
00204 inline G4LogicalVolume* BDSAcceleratorComponent::GetReadOutLogicalVolume() const
00205 {return readOutLV;}
00206
00207 inline G4double BDSAcceleratorComponent::GetParameterValue(G4String spec, G4String name) const
00208 {
00209 G4double value = 0;
00210
00211 std::string delimiters = "&";
00212 std::string param = name + "=";
00213
00214 int pos = spec.find(param);
00215 if( pos >= 0 )
00216 {
00217 int pos2 = spec.find("&",pos);
00218 int pos3 = spec.length();
00219 int tend = pos2 < 0 ? pos3 : pos2;
00220 int llen = tend - pos - param.length();
00221
00222 std::string val = spec.substr(pos + param.length(), llen);
00223
00224 value = atof(val.c_str());
00225 }
00226 return value;
00227 }
00228
00229 inline G4String BDSAcceleratorComponent::GetParameterValueString(G4String spec, G4String name) const
00230 {
00231 G4String value = "";
00232
00233 std::string delimiters = "&";
00234 std::string param = name + "=";
00235
00236 int pos = spec.find(param);
00237 if( pos >= 0 )
00238 {
00239 int pos2 = spec.find("&",pos);
00240 int pos3 = spec.length();
00241 int tend = pos2 < 0 ? pos3 : pos2;
00242 int llen = tend - pos - param.length();
00243
00244 value = spec.substr(pos + param.length(), llen);
00245 }
00246 return value;
00247 }
00248
00249 #endif