CaptDriftRegionBuilder.cc
Go to the documentation of this file.
3 #include "EDepSimBuilder.hh"
4 #include "EDepSimEMFieldSetup.hh"
5 #include "EDepSimUniformField.hh"
6 
7 #include <EDepSimLog.hh>
8 
9 #include <globals.hh>
10 #include <G4Material.hh>
11 #include <G4LogicalVolume.hh>
12 #include <G4VPhysicalVolume.hh>
13 #include <G4PVPlacement.hh>
14 #include <G4RotationMatrix.hh>
15 #include <G4VisAttributes.hh>
16 #include <G4UserLimits.hh>
17 #include <G4FieldManager.hh>
18 
19 #include <G4SystemOfUnits.hh>
20 #include <G4PhysicalConstants.hh>
21 
22 #include <G4Polyhedra.hh>
23 
25  : public EDepSim::BuilderMessenger {
26 private:
28  G4UIcmdWithADoubleAndUnit* fApothemCMD;
29  G4UIcmdWithADoubleAndUnit* fDriftLengthCMD;
30  G4UIcmdWithADoubleAndUnit* fWirePlaneSpacingCMD;
31 
32 public:
34  : EDepSim::BuilderMessenger(c,"Control the drift region geometry."),
35  fBuilder(c) {
36 
37  fApothemCMD
38  = new G4UIcmdWithADoubleAndUnit(CommandName("apothem"),this);
39  fApothemCMD->SetGuidance("Set the apothem of the drift region.");
40  fApothemCMD->SetParameterName("apothem",false);
41  fApothemCMD->SetUnitCategory("Length");
42 
43  fDriftLengthCMD = new G4UIcmdWithADoubleAndUnit(
44  CommandName("driftLength"),this);
45  fDriftLengthCMD->SetGuidance(
46  "Set the drift length from cathode to wires.");
47  fDriftLengthCMD->SetParameterName("drift",false);
48  fDriftLengthCMD->SetUnitCategory("Length");
49 
50  fWirePlaneSpacingCMD = new G4UIcmdWithADoubleAndUnit(
51  CommandName("wirePlaneSpacing"),this);
52  fWirePlaneSpacingCMD->SetGuidance(
53  "Set spacing between the wires.");
54  fWirePlaneSpacingCMD->SetParameterName("space",false);
55  fWirePlaneSpacingCMD->SetUnitCategory("Length");
56  };
57 
59  delete fApothemCMD;
60  delete fDriftLengthCMD;
61  delete fWirePlaneSpacingCMD;
62  };
63 
64  void SetNewValue(G4UIcommand *cmd, G4String val) {
65  if (cmd==fApothemCMD) {
66  fBuilder->SetApothem(fApothemCMD->GetNewDoubleValue(val));
67  }
68  else if (cmd==fDriftLengthCMD) {
69  fBuilder->SetDriftLength(fDriftLengthCMD->GetNewDoubleValue(val));
70  }
71  else if (cmd==fWirePlaneSpacingCMD) {
72  fBuilder->SetWirePlaneSpacing(
73  fWirePlaneSpacingCMD->GetNewDoubleValue(val));
74  }
75  else {
77  }
78  };
79 };
80 
82  SetMessenger(new CaptDriftRegionMessenger(this));
83  SetApothem(1037.5*mm);
84  SetDriftLength(1000*mm);
85  SetWirePlaneSpacing(3.18*mm);
86  SetSensitiveDetector("drift","segment");
87  SetMaximumHitLength(1.0*mm);
88  SetMaximumHitSagitta(0.25*mm);
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 }
96 
98 
100  CaptWirePlaneBuilder& wires = Get<CaptWirePlaneBuilder>("XPlane");
101  return GetDriftLength() + GetGridPlaneOffset() + wires.GetHeight();
102 }
103 
105  double zOffset = GetHeight()/2.0;
106  zOffset -= GetGridPlaneOffset();
107  return G4ThreeVector(0.0, 0.0, zOffset);
108 }
109 
111  CaptWirePlaneBuilder& wires = Get<CaptWirePlaneBuilder>("XPlane");
112  return 0.0*GetWirePlaneSpacing() + wires.GetHeight()/2;
113 }
114 
116  CaptWirePlaneBuilder& wires = Get<CaptWirePlaneBuilder>("XPlane");
117  return 1.0*GetWirePlaneSpacing() + wires.GetHeight()/2;
118 }
119 
121  CaptWirePlaneBuilder& wires = Get<CaptWirePlaneBuilder>("XPlane");
122  return 2.0*GetWirePlaneSpacing() + wires.GetHeight()/2;
123 }
124 
126  CaptWirePlaneBuilder& wires = Get<CaptWirePlaneBuilder>("XPlane");
127  return 3.0*GetWirePlaneSpacing() + wires.GetHeight()/2;
128 }
129 
131  CaptWirePlaneBuilder& wires = Get<CaptWirePlaneBuilder>("XPlane");
132  return 4.0*GetWirePlaneSpacing() + wires.GetHeight()/2;
133 }
134 
135 G4LogicalVolume *CaptDriftRegionBuilder::GetPiece(void) {
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
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.
CaptDriftRegionMessenger(CaptDriftRegionBuilder *c)
void SetNewValue(G4UIcommand *cmd, G4String val)
G4UIcmdWithADoubleAndUnit * fApothemCMD
double GetXPlaneOffset()
Get the offset of the X plane from the top of the drift region.
void SetWirePlaneSpacing(double v)
Set the spacing between the wire planes.
G4UIcmdWithADoubleAndUnit * fWirePlaneSpacingCMD
G4UIcmdWithADoubleAndUnit * fDriftLengthCMD
virtual G4LogicalVolume * GetPiece(void)
double GetVPlaneOffset()
Get the offset of the X plane from the top of the drift region.
Construct a module from components.
Definition: TG4HitSegment.h:10
BuilderMessenger(EDepSim::Builder *c, const char *guide=NULL)
virtual G4LogicalVolume * GetPiece(void)
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
CaptDriftRegionBuilder * fBuilder
void SetNewValue(G4UIcommand *cmd, G4String val)
static constexpr double mm
Definition: Units.h:65
static constexpr double degree
Definition: Units.h:161
list cmd
Definition: getreco.py:22
G4String CommandName(G4String cmd)
Build a command name with the directory prefix.