Classes | Public Member Functions | Private Types | Private Member Functions | Private Attributes | List of all members
CaptCryostatBuilder Class Reference

#include <CaptCryostatBuilder.hh>

Inheritance diagram for CaptCryostatBuilder:
EDepSim::Builder

Classes

class  Point
 

Public Member Functions

 CaptCryostatBuilder (G4String name, EDepSim::Builder *parent)
 
virtual ~CaptCryostatBuilder ()
 
virtual G4LogicalVolume * GetPiece (void)
 
G4ThreeVector GetOffset ()
 
G4ThreeVector GetTPCOffset ()
 
double GetArgonDepth () const
 
void SetArgonDepth (double v)
 
double GetTPCDepth () const
 
void SetTPCDepth (double v)
 
std::string GetVesselType () const
 
void SetVesselType (std::string v)
 
- 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 Types

typedef std::vector< PointShape
 

Private Member Functions

void Init (void)
 
void DefineCAPTAINVessel ()
 Fill the vessel definition with values for CAPTAIN. More...
 
void DefineMiniCAPTAINVessel ()
 Fill the vessel definition with values for miniCAPTAIN. More...
 

Private Attributes

G4LogicalVolume * fUllageVolume
 
G4LogicalVolume * fLiquidVolume
 
std::string fVesselType
 The name of the vessel to be built (currently CAPTAIN or mCAPTAIN). More...
 
double fArgonDepth
 The distance from the top of the flange to the liquid argon. More...
 
double fTPCDepth
 The distance from the top of the flange to the TPC. More...
 
Shape fInnerVessel
 The inner vessel polycone points. More...
 
Shape fOuterVessel
 The outer vessel polycone points. More...
 
Shape fVesselEnvelope
 The vessel envelop. 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 cryostat. This builds the vessel, as well as the argon gas and liquid inside. This builds the vessel directly and fills it with argon gas. The liquid argon is build by a separate constructor which will contain the TPC.

Bug:
The current class is an extremely simplistic model. I think the cryostat is multi-walled, but I don't have the design documents right now. The invariants for this class are the inner diameter and inner heights which define the cold volume. The outer diameter and height should be calculated.

Definition at line 19 of file CaptCryostatBuilder.hh.

Member Typedef Documentation

Definition at line 97 of file CaptCryostatBuilder.hh.

Constructor & Destructor Documentation

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

Definition at line 21 of file CaptCryostatBuilder.hh.

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

Definition at line 94 of file CaptCryostatBuilder.cc.

94 {}

Member Function Documentation

void CaptCryostatBuilder::DefineCAPTAINVessel ( )
private

Fill the vessel definition with values for CAPTAIN.

Definition at line 104 of file CaptCryostatBuilder.cc.

104  {
105  SetArgonDepth(24*25.4*CLHEP::mm); // The design spec for CAPTAIN (24 inch).
106  SetTPCDepth(GetArgonDepth()+25*CLHEP::mm); // I made this one up...
107 
108  fInnerVessel.clear();
109 #include "captainInnerVessel.hxx"
110 
111  fOuterVessel.clear();
112 #include "captainOuterVessel.hxx"
113 
114  for (Shape::reverse_iterator p = fOuterVessel.rbegin();
115  p != fOuterVessel.rend(); ++p) {
116  fVesselEnvelope.push_back(*p);
117  }
118  for (Shape::reverse_iterator p = fInnerVessel.rbegin();
119  p != fInnerVessel.rend(); ++p) {
120  if (p->fZ >= fVesselEnvelope.back().fZ) continue;
121  fVesselEnvelope.push_back(*p);
122  }
124 
125 }
Shape fInnerVessel
The inner vessel polycone points.
void SetArgonDepth(double v)
p
Definition: test.py:223
void SetTPCDepth(double v)
Shape fOuterVessel
The outer vessel polycone points.
static unsigned int reverse(QString &chars, unsigned char *level, unsigned int a, unsigned int b)
Definition: qstring.cpp:11649
Shape fVesselEnvelope
The vessel envelop.
double GetArgonDepth() const
static constexpr double mm
Definition: Units.h:65
void CaptCryostatBuilder::DefineMiniCAPTAINVessel ( )
private

Fill the vessel definition with values for miniCAPTAIN.

Definition at line 127 of file CaptCryostatBuilder.cc.

127  {
128  SetArgonDepth(9*25.4*CLHEP::mm); // The design spec for CAPTAIN (24 inch).
129  SetTPCDepth(GetArgonDepth()+25*CLHEP::mm); // I made this one up...
130 
131  fInnerVessel.clear();
133 
134  fOuterVessel.clear();
136 
137  for (Shape::reverse_iterator p = fOuterVessel.rbegin();
138  p != fOuterVessel.rend(); ++p) {
139  fVesselEnvelope.push_back(*p);
140  }
141  for (Shape::reverse_iterator p = fInnerVessel.rbegin();
142  p != fInnerVessel.rend(); ++p) {
143  if (p->fZ >= fVesselEnvelope.back().fZ) continue;
144  fVesselEnvelope.push_back(*p);
145  }
147 
148 }
Shape fInnerVessel
The inner vessel polycone points.
void SetArgonDepth(double v)
p
Definition: test.py:223
void SetTPCDepth(double v)
Shape fOuterVessel
The outer vessel polycone points.
static unsigned int reverse(QString &chars, unsigned char *level, unsigned int a, unsigned int b)
Definition: qstring.cpp:11649
Shape fVesselEnvelope
The vessel envelop.
double GetArgonDepth() const
static constexpr double mm
Definition: Units.h:65
double CaptCryostatBuilder::GetArgonDepth ( ) const
inline

Get (or set) the distance from the flange to the top of the liquid argon.

Definition at line 34 of file CaptCryostatBuilder.hh.

34 {return fArgonDepth;}
double fArgonDepth
The distance from the top of the flange to the liquid argon.
G4ThreeVector CaptCryostatBuilder::GetOffset ( )

Return the offset of the intended origin of the volume relative to the center of the logical volume. To get the origin at the right location (say originPosition), the logical volume should be positioned at originPosition-GetOffset(). The origin of the cryostat is chosen so that it's at the center of flange (the top of the "bucket" flange" defines the Z origin of the cryostat). The bottom of the wire plane assembly is at the origin (a decision will be made in the future as to whether this is the bottom of the grid, or the bottom of the V plane. This means that the wires for the V plane are at a (very) small positive z coordinate.

Definition at line 96 of file CaptCryostatBuilder.cc.

96  {
97  return G4ThreeVector(0,0,0);
98 }
G4LogicalVolume * CaptCryostatBuilder::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 150 of file CaptCryostatBuilder.cc.

150  {
151 
152  if (fVesselType == "CAPTAIN") {
154  }
155  else if (fVesselType == "mCAPTAIN") {
157  }
158  else {
159  std::cout << "Undefine vessel type: " << fVesselType << std::endl;
160  std::exit(0);
161  }
162 
163  std::vector<double> conePlanes;
164  std::vector<double> coneMax;
165  std::vector<double> coneMin;
166 
167  ////////////////////////////////////////////////////////
168  // Define the envelope to contain the vessel.
169  ////////////////////////////////////////////////////////
170  for (Shape::reverse_iterator p = fVesselEnvelope.rbegin();
171  p != fVesselEnvelope.rend();
172  ++p) {
173  conePlanes.push_back(- p->fZ);
174  coneMax.push_back(p->fOuter+10*CLHEP::cm);
175  coneMin.push_back(0.0);
176  }
177 
178  G4LogicalVolume* logVolume
179  = new G4LogicalVolume(
180  new G4Polycone(GetName(),
181  0*CLHEP::degree, 360*CLHEP::degree, conePlanes.size(),
182  conePlanes.data(),
183  coneMin.data(), coneMax.data()),
184  FindMaterial("Air"),
185  GetName());
186  logVolume->SetVisAttributes(G4VisAttributes::Invisible);
187 
188  ////////////////////////////////////////////////////////
189  // Define the outer vessel.
190  ////////////////////////////////////////////////////////
191  conePlanes.clear();
192  coneMax.clear();
193  coneMin.clear();
194  for (Shape::reverse_iterator p = fOuterVessel.rbegin();
195  p != fOuterVessel.rend();
196  ++p) {
197  conePlanes.push_back(- p->fZ);
198  coneMin.push_back(p->fInner);
199  coneMax.push_back(p->fOuter);
200  }
201 
202  G4LogicalVolume* logOuterVessel
203  = new G4LogicalVolume(
204  new G4Polycone(GetName()+"/OuterVessel",
205  0*CLHEP::degree, 360*CLHEP::degree, conePlanes.size(),
206  conePlanes.data(),
207  coneMin.data(), coneMax.data()),
208  FindMaterial("SS_304"),
209  GetName()+"/OuterVessel");
210  logOuterVessel->SetVisAttributes(GetColor(logOuterVessel));
211 
212  new G4PVPlacement(NULL, // rotation.
213  G4ThreeVector(0,0,0),
214  logOuterVessel, // logical volume
215  logOuterVessel->GetName(), // name
216  logVolume, // mother volume
217  false, // (not used)
218  0, // Copy number (zero)
219  Check()); // Check overlaps.
220 
221 
222  ////////////////////////////////////////////////////////
223  // Define the inner vessel.
224  ////////////////////////////////////////////////////////
225  conePlanes.clear();
226  coneMax.clear();
227  coneMin.clear();
228  for (Shape::reverse_iterator p = fInnerVessel.rbegin();
229  p != fInnerVessel.rend();
230  ++p) {
231  conePlanes.push_back(- p->fZ);
232  coneMin.push_back(p->fInner);
233  coneMax.push_back(p->fOuter);
234  }
235 
236  G4LogicalVolume* logInnerVessel
237  = new G4LogicalVolume(
238  new G4Polycone(GetName()+"/InnerVessel",
239  0*CLHEP::degree, 360*CLHEP::degree, conePlanes.size(),
240  conePlanes.data(),
241  coneMin.data(), coneMax.data()),
242  FindMaterial("SS_304"),
243  GetName()+"/InnerVessel");
244  logInnerVessel->SetVisAttributes(GetColor(logInnerVessel));
245 
246  new G4PVPlacement(NULL, // rotation.
247  G4ThreeVector(0,0,0),
248  logInnerVessel, // logical volume
249  logInnerVessel->GetName(), // name
250  logVolume, // mother volume
251  false, // (not used)
252  0, // Copy number (zero)
253  Check()); // Check overlaps.
254 
255  ////////////////////////////////////////////////////////
256  // Define the liquid volume.
257  ////////////////////////////////////////////////////////
258  conePlanes.clear();
259  coneMax.clear();
260  coneMin.clear();
261  for (Shape::reverse_iterator p = fInnerVessel.rbegin();
262  p != fInnerVessel.rend();
263  ++p) {
264  if (p->fZ < GetArgonDepth()) continue;
265  conePlanes.push_back(- p->fZ);
266  coneMin.push_back(0.0);
267  coneMax.push_back(p->fInner);
268  }
269  if (conePlanes.back() < -GetArgonDepth()) {
270  conePlanes.push_back(- GetArgonDepth());
271  coneMin.push_back(0.0);
272  coneMax.push_back(coneMax.back());
273  }
274 
276  = new G4LogicalVolume(
277  new G4Polycone(GetName()+"/Liquid",
278  0*CLHEP::degree, 360*CLHEP::degree, conePlanes.size(),
279  conePlanes.data(),
280  coneMin.data(), coneMax.data()),
281  FindMaterial("Argon_Liquid"),
282  GetName()+"/Liquid");
283  fLiquidVolume->SetVisAttributes(GetColor(fLiquidVolume));
284 
285  new G4PVPlacement(NULL, // rotation.
286  G4ThreeVector(0,0,0),
287  fLiquidVolume, // logical volume
288  fLiquidVolume->GetName()+"/Liquid", // name
289  logVolume, // mother volume
290  false, // (not used)
291  0, // Copy number (zero)
292  Check()); // Check overlaps.
293 
294  if (fVesselType == "CAPTAIN") {
295  CaptImmersedBuilder& immersed = Get<CaptImmersedBuilder>("Immersed");
296  G4LogicalVolume* logImmersed = immersed.GetPiece();
297  G4ThreeVector p = GetTPCOffset() - immersed.GetOffset();
298  new G4PVPlacement(NULL, // rotation.
299  p, // position
300  logImmersed, // logical volume
301  logImmersed->GetName(), // name
302  fLiquidVolume, // mother volume
303  false, // (not used)
304  0, // Copy number (zero)
305  Check()); // Check overlaps.
306  }
307  else if (fVesselType == "mCAPTAIN") {
308  MiniCaptImmersedBuilder& immersed
309  = Get<MiniCaptImmersedBuilder>("mImmersed");
310  G4LogicalVolume* logImmersed = immersed.GetPiece();
311  G4ThreeVector p = GetTPCOffset() - immersed.GetOffset();
312  new G4PVPlacement(NULL, // rotation.
313  p, // position
314  logImmersed, // logical volume
315  logImmersed->GetName(), // name
316  fLiquidVolume, // mother volume
317  false, // (not used)
318  0, // Copy number (zero)
319  Check()); // Check overlaps.
320  }
321  else {
322  std::cout << "Undefined immersed volume" << std::endl;
323  std::exit(0);
324  }
325 
326 
327  ////////////////////////////////////////////////////////
328  // Define the ullage volume.
329  ////////////////////////////////////////////////////////
330  conePlanes.clear();
331  coneMax.clear();
332  coneMin.clear();
333  for (Shape::reverse_iterator p = fInnerVessel.rbegin();
334  p != fInnerVessel.rend();
335  ++p) {
336  if (p->fZ >= GetArgonDepth()-1*CLHEP::mm) continue;
337  if (conePlanes.empty() ) {
338  conePlanes.push_back(-GetArgonDepth()+1*CLHEP::um);
339  coneMin.push_back(0.0);
340  coneMax.push_back(p->fInner);
341  }
342  conePlanes.push_back(- p->fZ);
343  coneMin.push_back(0.0);
344  coneMax.push_back(p->fInner);
345  }
346 
348  = new G4LogicalVolume(
349  new G4Polycone(GetName()+"/Ullage",
350  0*CLHEP::degree, 360*CLHEP::degree, conePlanes.size(),
351  conePlanes.data(),
352  coneMin.data(), coneMax.data()),
353  FindMaterial("Argon_Gas"),
354  GetName()+"/Ullage");
355  fUllageVolume->SetVisAttributes(GetColor(fUllageVolume));
356 
357  new G4PVPlacement(NULL, // rotation.
358  G4ThreeVector(0,0,0),
359  fUllageVolume, // logical volume
360  fUllageVolume->GetName(), // name
361  logVolume, // mother volume
362  false, // (not used)
363  0, // Copy number (zero)
364  Check()); // Check overlaps.
365 
366  if (fVesselType == "CAPTAIN") {
367  CaptExposedBuilder& exposed = Get<CaptExposedBuilder>("Exposed");
368  G4LogicalVolume* logExposed = exposed.GetPiece();
369  G4ThreeVector p(0.0,0.0,-GetArgonDepth());
370  p -= exposed.GetOffset();
371  new G4PVPlacement(NULL, // rotation.
372  p, // position
373  logExposed, // logical volume
374  logExposed->GetName(), // name
375  fUllageVolume, // mother volume
376  false, // (not used)
377  0, // Copy number (zero)
378  Check()); // Check overlaps.
379  }
380  else if (fVesselType == "mCAPTAIN") {
381  MiniCaptExposedBuilder& exposed
382  = Get<MiniCaptExposedBuilder>("mExposed");
383  G4LogicalVolume* logExposed = exposed.GetPiece();
384  if (logExposed) {
385  G4ThreeVector p(0.0,0.0,-GetArgonDepth());
386  p -= exposed.GetOffset();
387  new G4PVPlacement(NULL, // rotation.
388  p, // position
389  logExposed, // logical volume
390  logExposed->GetName(), // name
391  fUllageVolume, // mother volume
392  false, // (not used)
393  0, // Copy number (zero)
394  Check()); // Check overlaps.
395  }
396  }
397  else {
398  std::cout << "Undefined exposed volume" << std::endl;
399  std::exit(0);
400  }
401 
402  return logVolume;
403 }
static constexpr double cm
Definition: Units.h:68
G4String GetName(void)
Return the base name of the object that this builds.
void DefineMiniCAPTAINVessel()
Fill the vessel definition with values for miniCAPTAIN.
G4ThreeVector GetOffset()
G4LogicalVolume * fUllageVolume
G4Material * FindMaterial(G4String m)
Shape fInnerVessel
The inner vessel polycone points.
std::string fVesselType
The name of the vessel to be built (currently CAPTAIN or mCAPTAIN).
p
Definition: test.py:223
void DefineCAPTAINVessel()
Fill the vessel definition with values for CAPTAIN.
G4ThreeVector GetTPCOffset()
Shape fOuterVessel
The outer vessel polycone points.
virtual G4LogicalVolume * GetPiece(void)
G4VisAttributes GetColor(G4LogicalVolume *volume, double opacity=0.0)
Shape fVesselEnvelope
The vessel envelop.
virtual G4LogicalVolume * GetPiece(void)
double GetArgonDepth() const
static constexpr double mm
Definition: Units.h:65
G4ThreeVector GetOffset()
static constexpr double degree
Definition: Units.h:161
virtual G4LogicalVolume * GetPiece(void)
G4LogicalVolume * fLiquidVolume
QTextStream & endl(QTextStream &s)
virtual G4LogicalVolume * GetPiece(void)
double CaptCryostatBuilder::GetTPCDepth ( ) const
inline

Get (or set) the distance from the flange to the TPC center.

Definition at line 40 of file CaptCryostatBuilder.hh.

40 {return fTPCDepth;}
double fTPCDepth
The distance from the top of the flange to the TPC.
G4ThreeVector CaptCryostatBuilder::GetTPCOffset ( )

The position of the TPC origin (center of the grid plane defining the drift region) relative to the center of the Cryostat flange. The TPC origin should be arranged so that it will be at center of the global coordinate system.

Definition at line 100 of file CaptCryostatBuilder.cc.

100  {
101  return G4ThreeVector(0,0,-GetTPCDepth());
102 }
double GetTPCDepth() const
std::string CaptCryostatBuilder::GetVesselType ( ) const
inline

Get (or set) the type of vessel that should be built.

Definition at line 46 of file CaptCryostatBuilder.hh.

46 {return fVesselType;}
std::string fVesselType
The name of the vessel to be built (currently CAPTAIN or mCAPTAIN).
void CaptCryostatBuilder::Init ( void  )
private

Definition at line 81 of file CaptCryostatBuilder.cc.

81  {
83 
84  SetVesselType("CAPTAIN");
85 
86  SetSensitiveDetector("cryo","segment");
87 
88  AddBuilder(new CaptImmersedBuilder("Immersed",this));
89  AddBuilder(new CaptExposedBuilder("Exposed",this));
90  AddBuilder(new MiniCaptImmersedBuilder("mImmersed",this));
91  AddBuilder(new MiniCaptExposedBuilder("mExposed",this));
92 }
void SetMessenger(G4UImessenger *m)
Set the messenger for this constructor.
void AddBuilder(EDepSim::Builder *c)
virtual void SetSensitiveDetector(G4VSensitiveDetector *s)
Set the sensitive detector for this component.
void SetVesselType(std::string v)
void CaptCryostatBuilder::SetArgonDepth ( double  v)
inline

Definition at line 35 of file CaptCryostatBuilder.hh.

35 {fArgonDepth = v;}
double fArgonDepth
The distance from the top of the flange to the liquid argon.
void CaptCryostatBuilder::SetTPCDepth ( double  v)
inline

Definition at line 41 of file CaptCryostatBuilder.hh.

41 {fTPCDepth = v;}
double fTPCDepth
The distance from the top of the flange to the TPC.
void CaptCryostatBuilder::SetVesselType ( std::string  v)
inline

Definition at line 47 of file CaptCryostatBuilder.hh.

47 {fVesselType=v;}
std::string fVesselType
The name of the vessel to be built (currently CAPTAIN or mCAPTAIN).

Member Data Documentation

double CaptCryostatBuilder::fArgonDepth
private

The distance from the top of the flange to the liquid argon.

Definition at line 85 of file CaptCryostatBuilder.hh.

Shape CaptCryostatBuilder::fInnerVessel
private

The inner vessel polycone points.

Definition at line 100 of file CaptCryostatBuilder.hh.

G4LogicalVolume* CaptCryostatBuilder::fLiquidVolume
private

The volume for the liquid argon. This is set when the cryostat vessel is constructed and fills the liquid part of the inner void. It can be passed to the daughter volumes to define the shape of the liquid.

Definition at line 79 of file CaptCryostatBuilder.hh.

Shape CaptCryostatBuilder::fOuterVessel
private

The outer vessel polycone points.

Definition at line 102 of file CaptCryostatBuilder.hh.

double CaptCryostatBuilder::fTPCDepth
private

The distance from the top of the flange to the TPC.

Definition at line 88 of file CaptCryostatBuilder.hh.

G4LogicalVolume* CaptCryostatBuilder::fUllageVolume
private

The volume for the ullage. This is set when the cryostat vessel is constructed and fills the gas part of the inner void. It can be passed to daughter volumes to define the shape of the ullage.

Definition at line 74 of file CaptCryostatBuilder.hh.

Shape CaptCryostatBuilder::fVesselEnvelope
private

The vessel envelop.

Definition at line 104 of file CaptCryostatBuilder.hh.

std::string CaptCryostatBuilder::fVesselType
private

The name of the vessel to be built (currently CAPTAIN or mCAPTAIN).

Definition at line 82 of file CaptCryostatBuilder.hh.


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