Public Member Functions | Private Member Functions | Private Attributes | List of all members
CaptDriftRegionBuilder Class Reference

#include <CaptDriftRegionBuilder.hh>

Inheritance diagram for CaptDriftRegionBuilder:
EDepSim::Builder

Public Member Functions

 CaptDriftRegionBuilder (G4String name, EDepSim::Builder *parent)
 
virtual ~CaptDriftRegionBuilder ()
 
virtual G4LogicalVolume * GetPiece (void)
 
double GetHeight ()
 Get the total height of the drift region. More...
 
void SetWirePlaneSpacing (double v)
 Set the spacing between the wire planes. More...
 
double GetWirePlaneSpacing ()
 
G4ThreeVector GetOffset ()
 
double GetXPlaneOffset ()
 Get the offset of the X plane from the top of the drift region. More...
 
double GetVPlaneOffset ()
 Get the offset of the X plane from the top of the drift region. More...
 
double GetUPlaneOffset ()
 Get the offset of the X plane from the top of the drift region. More...
 
double GetGridPlaneOffset ()
 Get the offset of the grid from the top of the drift region. More...
 
double GetGroundPlaneOffset ()
 
void SetApothem (double v)
 
double GetApothem ()
 
void SetRadius (double v)
 
double GetRadius () const
 
void SetDriftLength (double v)
 
double GetDriftLength ()
 
- Public Member Functions inherited from EDepSim::Builder
 Builder (G4String n, EDepSim::UserDetectorConstruction *c)
 
 Builder (G4String n, EDepSim::Builder *parent)
 
virtual ~Builder ()
 
G4String GetName (void)
 Return the base name of the object that this builds. More...
 
G4String GetLocalName (void)
 Return the base name of the object that this builds. More...
 
void SetLocalName (const G4String &name)
 Set the base name of the logical volume that this builds. More...
 
void SetOpacity (double v)
 Set the relative opacity of the constructed object. More...
 
double GetOpacity (void) const
 Get the relative opacity of the constructed object. More...
 
void SetDaughterOpacity (double v)
 Set the relative opacity of the object daughters. More...
 
EDepSim::UserDetectorConstructionGetConstruction (void)
 Return the detector construction that is building this piece. More...
 
virtual void SetSensitiveDetector (G4VSensitiveDetector *s)
 Set the sensitive detector for this component. More...
 
virtual G4VSensitiveDetector * GetSensitiveDetector (void)
 Get the sensitive detector for this component. More...
 
virtual void SetSensitiveDetector (G4String name, G4String type)
 Set the sensitive detector for this component by name. More...
 
virtual void SetMaximumHitSagitta (double sagitta)
 
virtual void SetMaximumHitLength (double length)
 Set the maximum length of a single hit segment. More...
 
G4UImessenger * GetMessenger (void)
 Return the messenger for this constructor. More...
 
void SetMessenger (G4UImessenger *m)
 Set the messenger for this constructor. More...
 
void AddBuilder (EDepSim::Builder *c)
 
template<class T >
T & Get (G4String n)
 
template<class T >
T * Find (G4String n)
 
bool Check ()
 
void SetCheck (bool v)
 Set the check value. More...
 

Private Member Functions

void Init (void)
 

Private Attributes

double fApothem
 The radius of the circle that can fit inside. More...
 
double fDriftLength
 The distance from the cathode to the first wire plane. More...
 
double fWirePlaneSpacing
 The distance between the wire planes. More...
 

Additional Inherited Members

- Protected Member Functions inherited from EDepSim::Builder
G4Material * FindMaterial (G4String m)
 
G4VisAttributes GetColor (G4LogicalVolume *volume, double opacity=0.0)
 
G4VisAttributes GetColor (G4Material *volume, double opacity=0.0)
 

Detailed Description

Construct the drift region of the TPC. The wire planes are located at the top of the drift region, with the first wire for each plane located at the most negative X coordinate. The U plane runs from negative Y to positive Y. The V plane runs from positive Y to negative Y. The planes are oriented from top to bottom X, U, V (in otherwords, the first plane seen by the electrons is the V plane). The X plane is the collection plane.

Definition at line 14 of file CaptDriftRegionBuilder.hh.

Constructor & Destructor Documentation

CaptDriftRegionBuilder::CaptDriftRegionBuilder ( G4String  name,
EDepSim::Builder parent 
)
inline

Definition at line 16 of file CaptDriftRegionBuilder.hh.

17  : EDepSim::Builder(name,parent) {Init();};
static QCString name
Definition: declinfo.cpp:673
CaptDriftRegionBuilder::~CaptDriftRegionBuilder ( )
virtual

Definition at line 97 of file CaptDriftRegionBuilder.cc.

97 {}

Member Function Documentation

double CaptDriftRegionBuilder::GetApothem ( )
inline

Definition at line 33 of file CaptDriftRegionBuilder.hh.

33 {return fApothem;}
double fApothem
The radius of the circle that can fit inside.
double CaptDriftRegionBuilder::GetDriftLength ( )
inline

Definition at line 48 of file CaptDriftRegionBuilder.hh.

48 {return fDriftLength;}
double fDriftLength
The distance from the cathode to the first wire plane.
double CaptDriftRegionBuilder::GetGridPlaneOffset ( )

Get the offset of the grid from the top of the drift region.

Definition at line 130 of file CaptDriftRegionBuilder.cc.

130  {
131  CaptWirePlaneBuilder& wires = Get<CaptWirePlaneBuilder>("XPlane");
132  return 4.0*GetWirePlaneSpacing() + wires.GetHeight()/2;
133 }
double CaptDriftRegionBuilder::GetGroundPlaneOffset ( )

Get the offset of the ground from the top of the drift region. The distance between the grid and the ground defines the drift distance for the electrons.

Definition at line 110 of file CaptDriftRegionBuilder.cc.

110  {
111  CaptWirePlaneBuilder& wires = Get<CaptWirePlaneBuilder>("XPlane");
112  return 0.0*GetWirePlaneSpacing() + wires.GetHeight()/2;
113 }
double CaptDriftRegionBuilder::GetHeight ( )

Get the total height of the drift region.

Definition at line 99 of file CaptDriftRegionBuilder.cc.

99  {
100  CaptWirePlaneBuilder& wires = Get<CaptWirePlaneBuilder>("XPlane");
101  return GetDriftLength() + GetGridPlaneOffset() + wires.GetHeight();
102 }
double GetGridPlaneOffset()
Get the offset of the grid from the top of the drift region.
G4ThreeVector CaptDriftRegionBuilder::GetOffset ( )

Return the offset of the intended origin of the volume relative to the center of the logical volume. To get the intended origin at the right location (say originPosition), the logical volume should be positioned at originPosition-GetOffset(). For the drift region, the intended origin is at the "center" of the bottom of the wire planes is at the global origin.

Definition at line 104 of file CaptDriftRegionBuilder.cc.

104  {
105  double zOffset = GetHeight()/2.0;
106  zOffset -= GetGridPlaneOffset();
107  return G4ThreeVector(0.0, 0.0, zOffset);
108 }
double GetHeight()
Get the total height of the drift region.
double GetGridPlaneOffset()
Get the offset of the grid from the top of the drift region.
G4LogicalVolume * CaptDriftRegionBuilder::GetPiece ( void  )
virtual

Construct and return a G4 volume for the object. This is a pure virtual function, which means it must be implemented by the inheriting classes. This returns an unplaced logical volume which faces along the Z axis.

Implements EDepSim::Builder.

Definition at line 135 of file CaptDriftRegionBuilder.cc.

135  {
136 
137  double rInner[] = {0.0, 0.0};
138  double rOuter[] = {GetApothem(), GetApothem()};
139  double zPlane[] = {GetHeight()/2, -GetHeight()/2};
140 
141  EDepSimLog("Construct " << GetName()
142  << " with " << GetHeight()/mm << " mm height"
143  << " and " << GetDriftLength()/mm << " mm drift");
144 
145  G4LogicalVolume *logVolume
146  = new G4LogicalVolume(new G4Polyhedra(GetName(),
147  90*degree, 450*degree,
148  6, 2,
149  zPlane, rInner, rOuter),
150  FindMaterial("Argon_Liquid"),
151  GetName());
152  logVolume->SetVisAttributes(GetColor(logVolume));
153  G4FieldManager* manager = new G4FieldManager();
154  G4ElectroMagneticField* field = new EDepSim::UniformField(
155  G4ThreeVector(0,0,0), G4ThreeVector(0.0,0.0,500.0*volt/cm));
156  EDepSim::EMFieldSetup* fieldSetup
157  = new EDepSim::EMFieldSetup(field,manager);
158  if (!fieldSetup) {
159  EDepSimError("Field not properly created");
160  throw std::runtime_error("Field not properly created");
161  }
162  logVolume->SetFieldManager(manager,true);
163  if (GetSensitiveDetector()) {
164  logVolume->SetSensitiveDetector(GetSensitiveDetector());
165  logVolume->SetUserLimits(new G4UserLimits(1.0*mm));
166  }
167 
168  // The wire planes are rotated so that the local Z axis points into the
169  // drift volume, and the local X axis points along the original X axis
170  // (this puts the local Y axis along the -Y global axis.
171  G4RotationMatrix* xRotation = new G4RotationMatrix();
172  xRotation->rotateX(180*degree);
173 
174  CaptWirePlaneBuilder& xPlane = Get<CaptWirePlaneBuilder>("XPlane");
175  xPlane.SetApothem(GetApothem());
176  G4LogicalVolume *logX = xPlane.GetPiece();
177  new G4PVPlacement(xRotation, // rotation.
178  G4ThreeVector(0,0,
179  (GetHeight()/2 - GetXPlaneOffset())),
180  logX, // logical volume
181  logX->GetName(), // name
182  logVolume, // mother volume
183  false, // (not used)
184  0, // Copy number (zero)
185  Check()); // Check overlaps.
186 
187 
188  G4RotationMatrix* vRotation = new G4RotationMatrix();
189  vRotation->rotateX(180*degree);
190  vRotation->rotateZ(-60*degree);
191 
192  CaptWirePlaneBuilder& vPlane = Get<CaptWirePlaneBuilder>("VPlane");
193  vPlane.SetApothem(GetApothem());
194  G4LogicalVolume *logV = vPlane.GetPiece();
195  new G4PVPlacement(vRotation, // rotation.
196  G4ThreeVector(0,0,
197  (GetHeight()/2 - GetVPlaneOffset())),
198  logV, // logical volume
199  logV->GetName(), // name
200  logVolume, // mother volume
201  false, // (not used)
202  0, // Copy number (zero)
203  Check()); // Check overlaps.
204 
205 
206  G4RotationMatrix* uRotation = new G4RotationMatrix();
207  uRotation->rotateX(180*degree);
208  uRotation->rotateZ(60*degree);
209 
210  CaptWirePlaneBuilder& uPlane = Get<CaptWirePlaneBuilder>("UPlane");
211  uPlane.SetApothem(GetApothem());
212  G4LogicalVolume *logU = uPlane.GetPiece();
213  new G4PVPlacement(uRotation, // rotation.
214  G4ThreeVector(0,0,
215  (GetHeight()/2 - GetUPlaneOffset())),
216  logU, // logical volume
217  logU->GetName(), // name
218  logVolume, // mother volume
219  false, // (not used)
220  0, // Copy number (zero)
221  Check()); // Check overlaps.
222 
223  CaptWirePlaneBuilder& gridPlane = Get<CaptWirePlaneBuilder>("GridPlane");
224  gridPlane.SetApothem(GetApothem());
225  G4LogicalVolume *logGrid = gridPlane.GetPiece();
226  new G4PVPlacement(xRotation, // rotation.
227  G4ThreeVector(0,0,
228  (GetHeight()/2 - GetGridPlaneOffset())),
229  logGrid, // logical volume
230  logGrid->GetName(), // name
231  logVolume, // mother volume
232  false, // (not used)
233  0, // Copy number (zero)
234  Check()); // Check overlaps.
235 
236  CaptWirePlaneBuilder& grndPlane = Get<CaptWirePlaneBuilder>("GroundPlane");
237  grndPlane.SetApothem(GetApothem());
238  G4LogicalVolume *logGrnd = grndPlane.GetPiece();
239  new G4PVPlacement(xRotation, // rotation.
240  G4ThreeVector(0,0,
241  (GetHeight()/2 - GetGroundPlaneOffset())),
242  logGrnd, // logical volume
243  logGrnd->GetName(), // name
244  logVolume, // mother volume
245  false, // (not used)
246  0, // Copy number (zero)
247  Check()); // Check overlaps.
248 
249  return logVolume;
250 }
static constexpr double cm
Definition: Units.h:68
#define EDepSimLog(outStream)
Definition: EDepSimLog.hh:717
G4String GetName(void)
Return the base name of the object that this builds.
double GetHeight()
Get the total height of the drift region.
double GetGridPlaneOffset()
Get the offset of the grid from the top of the drift region.
virtual G4VSensitiveDetector * GetSensitiveDetector(void)
Get the sensitive detector for this component.
double GetXPlaneOffset()
Get the offset of the X plane from the top of the drift region.
G4Material * FindMaterial(G4String m)
virtual G4LogicalVolume * GetPiece(void)
double GetVPlaneOffset()
Get the offset of the X plane from the top of the drift region.
volt_as<> volt
Type of potential stored in volts, in double precision.
double GetUPlaneOffset()
Get the offset of the X plane from the top of the drift region.
#define EDepSimError(outStream)
Definition: EDepSimLog.hh:503
G4VisAttributes GetColor(G4LogicalVolume *volume, double opacity=0.0)
static constexpr double mm
Definition: Units.h:65
static constexpr double degree
Definition: Units.h:161
double CaptDriftRegionBuilder::GetRadius ( ) const
inline

Definition at line 41 of file CaptDriftRegionBuilder.hh.

41 {return fApothem/0.866025403784;}
double fApothem
The radius of the circle that can fit inside.
double CaptDriftRegionBuilder::GetUPlaneOffset ( )

Get the offset of the X plane from the top of the drift region.

Definition at line 125 of file CaptDriftRegionBuilder.cc.

125  {
126  CaptWirePlaneBuilder& wires = Get<CaptWirePlaneBuilder>("XPlane");
127  return 3.0*GetWirePlaneSpacing() + wires.GetHeight()/2;
128 }
double CaptDriftRegionBuilder::GetVPlaneOffset ( )

Get the offset of the X plane from the top of the drift region.

Definition at line 120 of file CaptDriftRegionBuilder.cc.

120  {
121  CaptWirePlaneBuilder& wires = Get<CaptWirePlaneBuilder>("XPlane");
122  return 2.0*GetWirePlaneSpacing() + wires.GetHeight()/2;
123 }
double CaptDriftRegionBuilder::GetWirePlaneSpacing ( )
inline

Definition at line 53 of file CaptDriftRegionBuilder.hh.

53 {return fWirePlaneSpacing;}
double fWirePlaneSpacing
The distance between the wire planes.
double CaptDriftRegionBuilder::GetXPlaneOffset ( )

Get the offset of the X plane from the top of the drift region.

Definition at line 115 of file CaptDriftRegionBuilder.cc.

115  {
116  CaptWirePlaneBuilder& wires = Get<CaptWirePlaneBuilder>("XPlane");
117  return 1.0*GetWirePlaneSpacing() + wires.GetHeight()/2;
118 }
void CaptDriftRegionBuilder::Init ( void  )
private

Definition at line 81 of file CaptDriftRegionBuilder.cc.

81  {
83  SetApothem(1037.5*mm);
84  SetDriftLength(1000*mm);
85  SetWirePlaneSpacing(3.18*mm);
86  SetSensitiveDetector("drift","segment");
89 
90  AddBuilder(new CaptWirePlaneBuilder("XPlane",this));
91  AddBuilder(new CaptWirePlaneBuilder("VPlane",this));
92  AddBuilder(new CaptWirePlaneBuilder("UPlane",this));
93  AddBuilder(new CaptWirePlaneBuilder("GridPlane",this));
94  AddBuilder(new CaptWirePlaneBuilder("GroundPlane",this));
95 }
virtual void SetMaximumHitLength(double length)
Set the maximum length of a single hit segment.
void SetMessenger(G4UImessenger *m)
Set the messenger for this constructor.
void SetWirePlaneSpacing(double v)
Set the spacing between the wire planes.
virtual void SetMaximumHitSagitta(double sagitta)
void AddBuilder(EDepSim::Builder *c)
virtual void SetSensitiveDetector(G4VSensitiveDetector *s)
Set the sensitive detector for this component.
static constexpr double mm
Definition: Units.h:65
void CaptDriftRegionBuilder::SetApothem ( double  v)
inline

Set the radius of the largest cylinder that can be contained in the drift region.

Definition at line 32 of file CaptDriftRegionBuilder.hh.

32 {fApothem = v;}
double fApothem
The radius of the circle that can fit inside.
void CaptDriftRegionBuilder::SetDriftLength ( double  v)
inline

Set the length of the drift region from the cathode to the first wire plane. (i.e. the grid plane).

Definition at line 47 of file CaptDriftRegionBuilder.hh.

47 {fDriftLength = v;}
double fDriftLength
The distance from the cathode to the first wire plane.
void CaptDriftRegionBuilder::SetRadius ( double  v)
inline

Set the radius of the smallest cylinder that contains the drift region. This is the maximum local Y dimension and is the apothem divided by the cosine of 30 degrees.

Definition at line 40 of file CaptDriftRegionBuilder.hh.

40 {fApothem = v*0.866025403784;}
double fApothem
The radius of the circle that can fit inside.
void CaptDriftRegionBuilder::SetWirePlaneSpacing ( double  v)
inline

Set the spacing between the wire planes.

Definition at line 52 of file CaptDriftRegionBuilder.hh.

52 {fWirePlaneSpacing = v;}
double fWirePlaneSpacing
The distance between the wire planes.

Member Data Documentation

double CaptDriftRegionBuilder::fApothem
private

The radius of the circle that can fit inside.

Definition at line 85 of file CaptDriftRegionBuilder.hh.

double CaptDriftRegionBuilder::fDriftLength
private

The distance from the cathode to the first wire plane.

Definition at line 88 of file CaptDriftRegionBuilder.hh.

double CaptDriftRegionBuilder::fWirePlaneSpacing
private

The distance between the wire planes.

Definition at line 91 of file CaptDriftRegionBuilder.hh.


The documentation for this class was generated from the following files: