CaptWirePlaneBuilder.cc
Go to the documentation of this file.
2 #include "EDepSimBuilder.hh"
3 
4 #include "EDepSimLog.hh"
5 
6 #include <globals.hh>
7 #include <G4Material.hh>
8 #include <G4LogicalVolume.hh>
9 #include <G4VPhysicalVolume.hh>
10 #include <G4PVPlacement.hh>
11 #include <G4VisAttributes.hh>
12 
13 #include <G4SystemOfUnits.hh>
14 #include <G4PhysicalConstants.hh>
15 
16 #include <G4Polyhedra.hh>
17 #include <G4Box.hh>
18 #include <G4Tubs.hh>
19 
20 #include <cmath>
21 
23  : public EDepSim::BuilderMessenger {
24 private:
26  G4UIcmdWithADoubleAndUnit* fApothemCMD;
27  G4UIcmdWithADoubleAndUnit* fSpacingCMD;
28  G4UIcmdWithAnInteger* fMaxWireCountCMD;
29 
30 
31 public:
33  : EDepSim::BuilderMessenger(c,"Control the drift region geometry."),
34  fBuilder(c) {
35 
36  fApothemCMD
37  = new G4UIcmdWithADoubleAndUnit(CommandName("apothem"),this);
38  fApothemCMD->SetGuidance("Set the apothem of the drift region.");
39  fApothemCMD->SetParameterName("apothem",false);
40  fApothemCMD->SetUnitCategory("Length");
41 
42  fSpacingCMD = new G4UIcmdWithADoubleAndUnit(
43  CommandName("wireSpacing"),this);
44  fSpacingCMD->SetGuidance(
45  "Set the wire spacing.");
46  fSpacingCMD->SetParameterName("spacing",false);
47  fSpacingCMD->SetUnitCategory("Length");
48 
49  fMaxWireCountCMD = new G4UIcmdWithAnInteger(
50  CommandName("maxWireCount"),this);
51  fMaxWireCountCMD->SetGuidance(
52  "Set the maximum number of wires in a plane.");
53  fMaxWireCountCMD->SetParameterName("count",false);
54 
55  }
56 
58  delete fApothemCMD;
59  delete fSpacingCMD;
60  delete fMaxWireCountCMD;
61  }
62 
63  void SetNewValue(G4UIcommand *cmd, G4String val) {
64  if (cmd==fApothemCMD) {
65  fBuilder->SetApothem(fApothemCMD->GetNewDoubleValue(val));
66  }
67  else if (cmd==fSpacingCMD) {
68  fBuilder->SetSpacing(fSpacingCMD->GetNewDoubleValue(val));
69  }
70  else if (cmd==fMaxWireCountCMD) {
71  fBuilder->SetMaxWireCount(fMaxWireCountCMD->GetNewIntValue(val));
72  }
73  else {
75  }
76  }
77 
78 };
79 
81  SetMessenger(new CaptWirePlaneMessenger(this));
82  SetApothem(1000*mm);
83  SetSpacing(3*mm);
84  SetHeight(1*mm);
85  SetMaxWireCount(0); // default value used to see if limit is set.
86  SetSensitiveDetector("drift","segment");
87  SetMaximumHitLength(1*mm);
88  SetMaximumHitSagitta(0.5*CLHEP::mm);
89 }
90 
92 
93 G4LogicalVolume *CaptWirePlaneBuilder::GetPiece(void) {
94 
95  double rInner[] = {0.0, 0.0};
96  double rOuter[] = {GetApothem(), GetApothem()};
97  double zPlane[] = {-GetHeight()/2, GetHeight()/2};
98 
99  G4LogicalVolume *logVolume
100  = new G4LogicalVolume(new G4Polyhedra(GetName(),
102  6, 2,
103  zPlane, rInner, rOuter),
104  FindMaterial("Argon_Liquid"),
105  GetName());
106  logVolume->SetVisAttributes(GetColor(FindMaterial("Captain_Wire"), 0.2));
107 
108  if (GetSensitiveDetector()) {
109  logVolume->SetSensitiveDetector(GetSensitiveDetector());
110  }
111 
112  // Determine the maximum (possible) wire length.
113  double maxLength = 2*GetApothem()/std::cos(30*CLHEP::degree);
114 
115  // Determine the number of wires. There is always at least one wire.
116  int wires = int (2.0*GetApothem()/GetSpacing());
117  if (wires<1) wires = 1;
118 
119  if (GetMaxWireCount()>0 && GetMaxWireCount()<wires) {
120  EDepSimLog("Reduce number of wires for " << GetName()
121  << " from " << wires
122  << " to " << GetMaxWireCount());
123  wires = GetMaxWireCount();
124  }
125 
126  EDepSimLog("Construct " << GetName()
127  << " with " << wires << " wires"
128  << " (spacing: " << GetSpacing()/CLHEP::mm << " mm)");
129 
130  if (wires<GetMaxWireCount()) {
131  EDepSimWarn(GetName() << ": Wire count is less than expected:"
132  << " " << wires
133  << " < " << GetMaxWireCount());
134  }
135 
136  G4RotationMatrix* wireRotation = NULL;
137 
138  // The core needs to be rotated into the wire volume. The core is not
139  // used inside the geometry.
140  G4RotationMatrix* coreRotation = new G4RotationMatrix();
141  coreRotation->rotateX(90*CLHEP::degree);
142 
143  // The offset of the first wire.
144  double baseOffset = - 0.5*GetSpacing()*(wires-1);
145  for (int wire = 0; wire<wires; ++wire) {
146  // The position of the wire.
147  double wireOffset = baseOffset + wire*GetSpacing();
148  double wireLength = maxLength-2*std::abs(wireOffset)*std::tan(30*CLHEP::degree);
149  double wireDiameter = std::min(GetHeight()/2, GetSpacing()/4);
150 
151  // Apply a correction to the length account for the corners of the
152  // wire box.
153  wireLength -= 2*GetSpacing()*std::cos(30*CLHEP::degree);
154 
155  G4LogicalVolume* logWire
156  = new G4LogicalVolume(new G4Box(GetName()+"/Wire",
157  GetSpacing()/2,
158  wireLength/2,
159  GetHeight()/2),
160  FindMaterial("Argon_Liquid"),
161  GetName()+"/Wire");
162  logWire->SetVisAttributes(G4VisAttributes::Invisible);
163 
164  if (GetSensitiveDetector()) {
165  logWire->SetSensitiveDetector(GetSensitiveDetector());
166  }
167 
168  new G4PVPlacement(wireRotation, // rotation.
169  G4ThreeVector(wireOffset,0,0), // position
170  logWire, // logical volume
171  logWire->GetName(), // name
172  logVolume, // mother volume
173  false, // (not used)
174  0, // Copy number (zero)
175  false); // Check overlaps.
176 
177  G4LogicalVolume* logCore
178  = new G4LogicalVolume(new G4Tubs(GetName()+"/Wire/Core",
179  0.0, wireDiameter/2.0,
180  wireLength/2,
182  FindMaterial("Captain_Wire"),
183  GetName()+"/Wire/Core");
184  // Draw this as if it was a wire plane.
185  logCore->SetVisAttributes(GetColor(logCore));
186 
187  new G4PVPlacement(coreRotation, // rotation.
188  G4ThreeVector(0,0,0), // position
189  logCore, // logical volume
190  logCore->GetName(), // name
191  logWire, // mother volume
192  false, // (not used)
193  0, // Copy number (zero)
194  false); // Check overlaps.
195 
196 
197 
198  }
199 
200  return logVolume;
201 }
#define EDepSimLog(outStream)
Definition: EDepSimLog.hh:717
void SetNewValue(G4UIcommand *cmd, G4String val)
void SetNewValue(G4UIcommand *cmd, G4String val)
CaptWirePlaneBuilder * fBuilder
G4UIcmdWithADoubleAndUnit * fSpacingCMD
virtual G4LogicalVolume * GetPiece(void)
T abs(T value)
Construct a module from components.
Definition: TG4HitSegment.h:10
CaptWirePlaneMessenger(CaptWirePlaneBuilder *c)
BuilderMessenger(EDepSim::Builder *c, const char *guide=NULL)
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55
#define EDepSimWarn(outStream)
Definition: EDepSimLog.hh:576
G4UIcmdWithAnInteger * fMaxWireCountCMD
G4UIcmdWithADoubleAndUnit * fApothemCMD
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.