Geant4  10.02.p03
LXeDetectorConstruction Class Reference

#include <LXeDetectorConstruction.hh>

Inheritance diagram for LXeDetectorConstruction:
Collaboration diagram for LXeDetectorConstruction:

Public Member Functions

 LXeDetectorConstruction ()
 
virtual ~LXeDetectorConstruction ()
 
virtual G4VPhysicalVolumeConstruct ()
 
virtual void ConstructSDandField ()
 
void SetDimensions (G4ThreeVector)
 
void SetHousingThickness (G4double)
 
void SetNX (G4int)
 
void SetNY (G4int)
 
void SetNZ (G4int)
 
void SetPMTRadius (G4double)
 
void SetDefaults ()
 
G4int GetNX ()
 
G4int GetNY ()
 
G4int GetNZ ()
 
G4double GetScintX ()
 
G4double GetScintY ()
 
G4double GetScintZ ()
 
G4double GetHousingThickness ()
 
G4double GetPMTRadius ()
 
G4double GetSlabZ ()
 
void SetSphereOn (G4bool)
 
void SetHousingReflectivity (G4double)
 
G4double GetHousingReflectivity ()
 
void SetWLSSlabOn (G4bool b)
 
G4bool GetWLSSlabOn ()
 
void SetMainVolumeOn (G4bool b)
 
G4bool GetMainVolumeOn ()
 
void SetNFibers (G4int n)
 
G4int GetNFibers ()
 
void SetMainScintYield (G4double)
 
void SetWLSScintYield (G4double)
 
- Public Member Functions inherited from G4VUserDetectorConstruction
 G4VUserDetectorConstruction ()
 
virtual ~G4VUserDetectorConstruction ()
 
virtual void CloneSD ()
 
virtual void CloneF ()
 
void RegisterParallelWorld (G4VUserParallelWorld *)
 
G4int ConstructParallelGeometries ()
 
void ConstructParallelSD ()
 
G4int GetNumberOfParallelWorld () const
 
G4VUserParallelWorldGetParallelWorld (G4int i) const
 

Static Public Member Functions

static G4bool GetSphereOn ()
 

Private Member Functions

void DefineMaterials ()
 
G4VPhysicalVolumeConstructDetector ()
 

Private Attributes

LXeDetectorMessengerfDetectorMessenger
 
G4BoxfExperimentalHall_box
 
G4LogicalVolumefExperimentalHall_log
 
G4VPhysicalVolumefExperimentalHall_phys
 
G4MaterialfLXe
 
G4MaterialfAl
 
G4ElementfN
 
G4ElementfO
 
G4MaterialfAir
 
G4MaterialfVacuum
 
G4ElementfC
 
G4ElementfH
 
G4MaterialfGlass
 
G4MaterialfPstyrene
 
G4MaterialfPMMA
 
G4MaterialfPethylene1
 
G4MaterialfPethylene2
 
G4double fScint_x
 
G4double fScint_y
 
G4double fScint_z
 
G4double fD_mtl
 
G4int fNx
 
G4int fNy
 
G4int fNz
 
G4double fOuterRadius_pmt
 
G4int fNfibers
 
G4double fRefl
 
G4bool fWLSslab
 
G4bool fMainVolumeOn
 
G4double fSlab_z
 
LXeMainVolumefMainVolume
 
G4MaterialPropertiesTablefLXe_mt
 
G4MaterialPropertiesTablefMPTPStyrene
 
G4Cache< LXeScintSD * > fScint_SD
 
G4Cache< LXePMTSD * > fPmt_SD
 

Static Private Attributes

static G4bool fSphereOn = true
 

Additional Inherited Members

- Protected Member Functions inherited from G4VUserDetectorConstruction
void SetSensitiveDetector (const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
 
void SetSensitiveDetector (G4LogicalVolume *logVol, G4VSensitiveDetector *aSD)
 

Detailed Description

Definition at line 53 of file LXeDetectorConstruction.hh.

Constructor & Destructor Documentation

◆ LXeDetectorConstruction()

LXeDetectorConstruction::LXeDetectorConstruction ( )

Definition at line 68 of file LXeDetectorConstruction.cc.

69 : fLXe_mt(NULL), fMPTPStyrene(NULL)
70 {
71  fExperimentalHall_box = NULL;
72  fExperimentalHall_log = NULL;
74 
75  fLXe = fAl = fAir = fVacuum = fGlass = NULL;
77 
78  fN = fO = fC = fH = NULL;
79 
80  SetDefaults();
81 
83 }
G4MaterialPropertiesTable * fMPTPStyrene
LXeDetectorMessenger * fDetectorMessenger
G4MaterialPropertiesTable * fLXe_mt
G4LogicalVolume * fExperimentalHall_log
G4VPhysicalVolume * fExperimentalHall_phys
Here is the call graph for this function:

◆ ~LXeDetectorConstruction()

LXeDetectorConstruction::~LXeDetectorConstruction ( )
virtual

Definition at line 87 of file LXeDetectorConstruction.cc.

87 {}

Member Function Documentation

◆ Construct()

G4VPhysicalVolume * LXeDetectorConstruction::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Definition at line 242 of file LXeDetectorConstruction.cc.

242  {
243 
251  }
252 
253  DefineMaterials();
254  return ConstructDetector();
255 }
static void Clean()
Definition: G4SolidStore.cc:79
static G4PhysicalVolumeStore * GetInstance()
G4VPhysicalVolume * ConstructDetector()
static G4LogicalVolumeStore * GetInstance()
static G4SolidStore * GetInstance()
static G4GeometryManager * GetInstance()
void OpenGeometry(G4VPhysicalVolume *vol=0)
G4VPhysicalVolume * fExperimentalHall_phys
Here is the call graph for this function:

◆ ConstructDetector()

G4VPhysicalVolume * LXeDetectorConstruction::ConstructDetector ( )
private

Definition at line 259 of file LXeDetectorConstruction.cc.

260 {
261  //The experimental hall walls are all 1m away from housing walls
265 
266  //Create experimental hall
268  = new G4Box("expHall_box",expHall_x,expHall_y,expHall_z);
270  fVacuum,"expHall_log",0,0,0);
272  fExperimentalHall_log,"expHall",0,false,0);
273 
275 
276  //Place the main volume
277  if(fMainVolumeOn){
279  = new LXeMainVolume(0,G4ThreeVector(),fExperimentalHall_log,false,0,this);
280  }
281 
282  //Place the WLS slab
283  if(fWLSslab){
284  G4VPhysicalVolume* slab = new LXeWLSSlab(0,G4ThreeVector(0.,0.,
285  -fScint_z/2.-fSlab_z-1.*cm),
286  fExperimentalHall_log,false,0,
287  this);
288 
289  //Surface properties for the WLS slab
290  G4OpticalSurface* scintWrap = new G4OpticalSurface("ScintWrap");
291 
292  new G4LogicalBorderSurface("ScintWrap", slab,
294  scintWrap);
295 
296  scintWrap->SetType(dielectric_metal);
297  scintWrap->SetFinish(polished);
298  scintWrap->SetModel(glisur);
299 
300  G4double pp[] = {2.0*eV, 3.5*eV};
301  const G4int num = sizeof(pp)/sizeof(G4double);
302  G4double reflectivity[] = {1., 1.};
303  assert(sizeof(reflectivity) == sizeof(pp));
304  G4double efficiency[] = {0.0, 0.0};
305  assert(sizeof(efficiency) == sizeof(pp));
306 
307  G4MaterialPropertiesTable* scintWrapProperty
309 
310  scintWrapProperty->AddProperty("REFLECTIVITY",pp,reflectivity,num);
311  scintWrapProperty->AddProperty("EFFICIENCY",pp,efficiency,num);
312  scintWrap->SetMaterialPropertiesTable(scintWrapProperty);
313  }
314 
315  return fExperimentalHall_phys;
316 }
void SetFinish(const G4OpticalSurfaceFinish)
static const double cm
Definition: G4SIunits.hh:118
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:64
int G4int
Definition: G4Types.hh:78
G4MaterialPropertyVector * AddProperty(const char *key, G4double *PhotonEnergies, G4double *PropertyValues, G4int NumEntries)
static const double eV
Definition: G4SIunits.hh:212
static const G4VisAttributes Invisible
G4LogicalVolume * fExperimentalHall_log
static const double m
Definition: G4SIunits.hh:128
double G4double
Definition: G4Types.hh:76
void SetModel(const G4OpticalSurfaceModel model)
void SetMaterialPropertiesTable(G4MaterialPropertiesTable *anMPT)
void SetType(const G4SurfaceType &type)
void SetVisAttributes(const G4VisAttributes *pVA)
G4VPhysicalVolume * fExperimentalHall_phys
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructSDandField()

void LXeDetectorConstruction::ConstructSDandField ( )
virtual

Reimplemented from G4VUserDetectorConstruction.

Definition at line 320 of file LXeDetectorConstruction.cc.

320  {
321 
322  if (!fMainVolume) return;
323 
324  // PMT SD
325 
326  if (!fPmt_SD.Get()) {
327  //Created here so it exists as pmts are being placed
328  G4cout << "Construction /LXeDet/pmtSD" << G4endl;
329  LXePMTSD* pmt_SD = new LXePMTSD("/LXeDet/pmtSD");
330  fPmt_SD.Put(pmt_SD);
331 
332  pmt_SD->InitPMTs((fNx*fNy+fNx*fNz+fNy*fNz)*2); //let pmtSD know # of pmts
334  }
335 
336  //sensitive detector is not actually on the photocathode.
337  //processHits gets done manually by the stepping action.
338  //It is used to detect when photons hit and get absorbed&detected at the
339  //boundary to the photocathode (which doesnt get done by attaching it to a
340  //logical volume.
341  //It does however need to be attached to something or else it doesnt get
342  //reset at the begining of events
343 
345 
346  // Scint SD
347 
348  if (!fScint_SD.Get()) {
349  G4cout << "Construction /LXeDet/scintSD" << G4endl;
350  LXeScintSD* scint_SD = new LXeScintSD("/LXeDet/scintSD");
351  fScint_SD.Put(scint_SD);
352  }
354 }
void Put(const value_type &val) const
Definition: G4Cache.hh:286
value_type & Get() const
Definition: G4Cache.hh:282
std::vector< G4ThreeVector > GetPmtPositions()
G4LogicalVolume * GetLogPhotoCath()
G4LogicalVolume * GetLogScint()
G4GLOB_DLL std::ostream G4cout
void SetSensitiveDetector(const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
void SetPmtPositions(const std::vector< G4ThreeVector > &positions)
Definition: LXePMTSD.cc:62
G4Cache< LXeScintSD * > fScint_SD
void InitPMTs(G4int nPMTs)
Definition: LXePMTSD.hh:64
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:

◆ DefineMaterials()

void LXeDetectorConstruction::DefineMaterials ( void  )
private

Definition at line 91 of file LXeDetectorConstruction.cc.

91  {
92  G4double a; // atomic mass
93  G4double z; // atomic number
95 
96  G4int polyPMMA = 1;
97  G4int nC_PMMA = 3+2*polyPMMA;
98  G4int nH_PMMA = 6+2*polyPMMA;
99 
100  G4int polyeth = 1;
101  G4int nC_eth = 2*polyeth;
102  G4int nH_eth = 4*polyeth;
103 
104  //***Elements
105  fH = new G4Element("H", "H", z=1., a=1.01*g/mole);
106  fC = new G4Element("C", "C", z=6., a=12.01*g/mole);
107  fN = new G4Element("N", "N", z=7., a= 14.01*g/mole);
108  fO = new G4Element("O" , "O", z=8., a= 16.00*g/mole);
109 
110  //***Materials
111  //Liquid Xenon
112  fLXe = new G4Material("LXe",z=54.,a=131.29*g/mole,density=3.020*g/cm3);
113  //Aluminum
114  fAl = new G4Material("Al",z=13.,a=26.98*g/mole,density=2.7*g/cm3);
115  //Vacuum
116  fVacuum = new G4Material("Vacuum",z=1.,a=1.01*g/mole,
118  1.e-19*pascal);
119  //Air
120  fAir = new G4Material("Air", density= 1.29*mg/cm3, 2);
121  fAir->AddElement(fN, 70*perCent);
122  fAir->AddElement(fO, 30*perCent);
123  //Glass
124  fGlass = new G4Material("Glass", density=1.032*g/cm3,2);
125  fGlass->AddElement(fC,91.533*perCent);
126  fGlass->AddElement(fH,8.467*perCent);
127  //Polystyrene
128  fPstyrene = new G4Material("Polystyrene", density= 1.03*g/cm3, 2);
129  fPstyrene->AddElement(fC, 8);
130  fPstyrene->AddElement(fH, 8);
131  //Fiber(PMMA)
132  fPMMA = new G4Material("PMMA", density=1190*kg/m3,3);
133  fPMMA->AddElement(fH,nH_PMMA);
134  fPMMA->AddElement(fC,nC_PMMA);
135  fPMMA->AddElement(fO,2);
136  //Cladding(polyethylene)
137  fPethylene1 = new G4Material("Pethylene1", density=1200*kg/m3,2);
138  fPethylene1->AddElement(fH,nH_eth);
139  fPethylene1->AddElement(fC,nC_eth);
140  //Double cladding(flourinated polyethylene)
141  fPethylene2 = new G4Material("Pethylene2", density=1400*kg/m3,2);
142  fPethylene2->AddElement(fH,nH_eth);
143  fPethylene2->AddElement(fC,nC_eth);
144 
145  //***Material properties tables
146 
147  G4double lxe_Energy[] = { 7.0*eV , 7.07*eV, 7.14*eV };
148  const G4int lxenum = sizeof(lxe_Energy)/sizeof(G4double);
149 
150  G4double lxe_SCINT[] = { 0.1, 1.0, 0.1 };
151  assert(sizeof(lxe_SCINT) == sizeof(lxe_Energy));
152  G4double lxe_RIND[] = { 1.59 , 1.57, 1.54 };
153  assert(sizeof(lxe_RIND) == sizeof(lxe_Energy));
154  G4double lxe_ABSL[] = { 35.*cm, 35.*cm, 35.*cm};
155  assert(sizeof(lxe_ABSL) == sizeof(lxe_Energy));
157  fLXe_mt->AddProperty("FASTCOMPONENT", lxe_Energy, lxe_SCINT, lxenum);
158  fLXe_mt->AddProperty("SLOWCOMPONENT", lxe_Energy, lxe_SCINT, lxenum);
159  fLXe_mt->AddProperty("RINDEX", lxe_Energy, lxe_RIND, lxenum);
160  fLXe_mt->AddProperty("ABSLENGTH", lxe_Energy, lxe_ABSL, lxenum);
161  fLXe_mt->AddConstProperty("SCINTILLATIONYIELD",12000./MeV);
162  fLXe_mt->AddConstProperty("RESOLUTIONSCALE",1.0);
163  fLXe_mt->AddConstProperty("FASTTIMECONSTANT",20.*ns);
164  fLXe_mt->AddConstProperty("SLOWTIMECONSTANT",45.*ns);
165  fLXe_mt->AddConstProperty("YIELDRATIO",1.0);
167 
168  // Set the Birks Constant for the LXe scintillator
169 
171 
172  G4double glass_RIND[]={1.49,1.49,1.49};
173  assert(sizeof(glass_RIND) == sizeof(lxe_Energy));
174  G4double glass_AbsLength[]={420.*cm,420.*cm,420.*cm};
175  assert(sizeof(glass_AbsLength) == sizeof(lxe_Energy));
177  glass_mt->AddProperty("ABSLENGTH",lxe_Energy,glass_AbsLength,lxenum);
178  glass_mt->AddProperty("RINDEX",lxe_Energy,glass_RIND,lxenum);
180 
181  G4double vacuum_Energy[]={2.0*eV,7.0*eV,7.14*eV};
182  const G4int vacnum = sizeof(vacuum_Energy)/sizeof(G4double);
183  G4double vacuum_RIND[]={1.,1.,1.};
184  assert(sizeof(vacuum_RIND) == sizeof(vacuum_Energy));
186  vacuum_mt->AddProperty("RINDEX", vacuum_Energy, vacuum_RIND,vacnum);
188  fAir->SetMaterialPropertiesTable(vacuum_mt);//Give air the same rindex
189 
190  G4double wls_Energy[] = {2.00*eV,2.87*eV,2.90*eV,3.47*eV};
191  const G4int wlsnum = sizeof(wls_Energy)/sizeof(G4double);
192 
193  G4double rIndexPstyrene[]={ 1.5, 1.5, 1.5, 1.5};
194  assert(sizeof(rIndexPstyrene) == sizeof(wls_Energy));
195  G4double absorption1[]={2.*cm, 2.*cm, 2.*cm, 2.*cm};
196  assert(sizeof(absorption1) == sizeof(wls_Energy));
197  G4double scintilFast[]={0.00, 0.00, 1.00, 1.00};
198  assert(sizeof(scintilFast) == sizeof(wls_Energy));
200  fMPTPStyrene->AddProperty("RINDEX",wls_Energy,rIndexPstyrene,wlsnum);
201  fMPTPStyrene->AddProperty("ABSLENGTH",wls_Energy,absorption1,wlsnum);
202  fMPTPStyrene->AddProperty("FASTCOMPONENT",wls_Energy, scintilFast,wlsnum);
203  fMPTPStyrene->AddConstProperty("SCINTILLATIONYIELD",10./keV);
204  fMPTPStyrene->AddConstProperty("RESOLUTIONSCALE",1.0);
205  fMPTPStyrene->AddConstProperty("FASTTIMECONSTANT", 10.*ns);
207 
208  // Set the Birks Constant for the Polystyrene scintillator
209 
211 
212  G4double RefractiveIndexFiber[]={ 1.60, 1.60, 1.60, 1.60};
213  assert(sizeof(RefractiveIndexFiber) == sizeof(wls_Energy));
214  G4double AbsFiber[]={9.00*m,9.00*m,0.1*mm,0.1*mm};
215  assert(sizeof(AbsFiber) == sizeof(wls_Energy));
216  G4double EmissionFib[]={1.0, 1.0, 0.0, 0.0};
217  assert(sizeof(EmissionFib) == sizeof(wls_Energy));
219  fiberProperty->AddProperty("RINDEX",wls_Energy,RefractiveIndexFiber,wlsnum);
220  fiberProperty->AddProperty("WLSABSLENGTH",wls_Energy,AbsFiber,wlsnum);
221  fiberProperty->AddProperty("WLSCOMPONENT",wls_Energy,EmissionFib,wlsnum);
222  fiberProperty->AddConstProperty("WLSTIMECONSTANT", 0.5*ns);
223  fPMMA->SetMaterialPropertiesTable(fiberProperty);
224 
225  G4double RefractiveIndexClad1[]={ 1.49, 1.49, 1.49, 1.49};
226  assert(sizeof(RefractiveIndexClad1) == sizeof(wls_Energy));
228  clad1Property->AddProperty("RINDEX",wls_Energy,RefractiveIndexClad1,wlsnum);
229  clad1Property->AddProperty("ABSLENGTH",wls_Energy,AbsFiber,wlsnum);
230  fPethylene1->SetMaterialPropertiesTable(clad1Property);
231 
232  G4double RefractiveIndexClad2[]={ 1.42, 1.42, 1.42, 1.42};
233  assert(sizeof(RefractiveIndexClad2) == sizeof(wls_Energy));
235  clad2Property->AddProperty("RINDEX",wls_Energy,RefractiveIndexClad2,wlsnum);
236  clad2Property->AddProperty("ABSLENGTH",wls_Energy,AbsFiber,wlsnum);
237  fPethylene2->SetMaterialPropertiesTable(clad2Property);
238 }
static const double cm
Definition: G4SIunits.hh:118
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
static const double MeV
Definition: G4SIunits.hh:211
static const double m3
Definition: G4SIunits.hh:130
int universe_mean_density
Definition: hepunit.py:307
void SetMaterialPropertiesTable(G4MaterialPropertiesTable *anMPT)
Definition: G4Material.hh:249
void SetBirksConstant(G4double value)
static const double mg
Definition: G4SIunits.hh:181
G4MaterialPropertiesTable * fMPTPStyrene
int G4int
Definition: G4Types.hh:78
G4MaterialPropertyVector * AddProperty(const char *key, G4double *PhotonEnergies, G4double *PropertyValues, G4int NumEntries)
G4double density
Definition: TRTMaterials.hh:39
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5206
#define pascal
static const double cm3
Definition: G4SIunits.hh:120
static const double kg
Definition: G4SIunits.hh:179
static const double perCent
Definition: G4SIunits.hh:329
void AddConstProperty(const char *key, G4double PropertyValue)
static const double kelvin
Definition: G4SIunits.hh:278
static const double eV
Definition: G4SIunits.hh:212
G4MaterialPropertiesTable * fLXe_mt
static const double mole
Definition: G4SIunits.hh:283
static const double m
Definition: G4SIunits.hh:128
void AddElement(G4Element *element, G4int nAtoms)
Definition: G4Material.cc:364
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
static const double mm
Definition: G4SIunits.hh:114
#define ns
Definition: xmlparse.cc:614
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetHousingReflectivity()

G4double LXeDetectorConstruction::GetHousingReflectivity ( )
inline

Definition at line 87 of file LXeDetectorConstruction.hh.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetHousingThickness()

G4double LXeDetectorConstruction::GetHousingThickness ( )
inline

Definition at line 79 of file LXeDetectorConstruction.hh.

Here is the caller graph for this function:

◆ GetMainVolumeOn()

G4bool LXeDetectorConstruction::GetMainVolumeOn ( )
inline

Definition at line 93 of file LXeDetectorConstruction.hh.

Here is the call graph for this function:

◆ GetNFibers()

G4int LXeDetectorConstruction::GetNFibers ( )
inline

Definition at line 96 of file LXeDetectorConstruction.hh.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetNX()

G4int LXeDetectorConstruction::GetNX ( )
inline

Definition at line 73 of file LXeDetectorConstruction.hh.

Here is the caller graph for this function:

◆ GetNY()

G4int LXeDetectorConstruction::GetNY ( )
inline

Definition at line 74 of file LXeDetectorConstruction.hh.

Here is the caller graph for this function:

◆ GetNZ()

G4int LXeDetectorConstruction::GetNZ ( )
inline

Definition at line 75 of file LXeDetectorConstruction.hh.

Here is the caller graph for this function:

◆ GetPMTRadius()

G4double LXeDetectorConstruction::GetPMTRadius ( )
inline

Definition at line 80 of file LXeDetectorConstruction.hh.

Here is the caller graph for this function:

◆ GetScintX()

G4double LXeDetectorConstruction::GetScintX ( )
inline

Definition at line 76 of file LXeDetectorConstruction.hh.

Here is the caller graph for this function:

◆ GetScintY()

G4double LXeDetectorConstruction::GetScintY ( )
inline

Definition at line 77 of file LXeDetectorConstruction.hh.

Here is the caller graph for this function:

◆ GetScintZ()

G4double LXeDetectorConstruction::GetScintZ ( )
inline

Definition at line 78 of file LXeDetectorConstruction.hh.

Here is the caller graph for this function:

◆ GetSlabZ()

G4double LXeDetectorConstruction::GetSlabZ ( )
inline

Definition at line 81 of file LXeDetectorConstruction.hh.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetSphereOn()

static G4bool LXeDetectorConstruction::GetSphereOn ( )
inlinestatic

Definition at line 84 of file LXeDetectorConstruction.hh.

84 {return fSphereOn;}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetWLSSlabOn()

G4bool LXeDetectorConstruction::GetWLSSlabOn ( )
inline

Definition at line 90 of file LXeDetectorConstruction.hh.

Here is the call graph for this function:

◆ SetDefaults()

void LXeDetectorConstruction::SetDefaults ( )

Definition at line 402 of file LXeDetectorConstruction.cc.

402  {
403 
404  //Resets to default values
405  fD_mtl=0.0635*cm;
406 
407  fScint_x = 17.8*cm;
408  fScint_y = 17.8*cm;
409  fScint_z = 22.6*cm;
410 
411  fNx = 2;
412  fNy = 2;
413  fNz = 3;
414 
415  fOuterRadius_pmt = 2.3*cm;
416 
417  fSphereOn = true;
418  fRefl=1.0;
419 
420  fNfibers=15;
421  fWLSslab=false;
422  fMainVolumeOn=true;
423  fMainVolume=NULL;
424  fSlab_z=2.5*mm;
425 
427  ->ApplyCommand("/LXe/detector/scintYieldFactor 1.");
428 
429  if(fLXe_mt)fLXe_mt->AddConstProperty("SCINTILLATIONYIELD",12000./MeV);
430  if(fMPTPStyrene)fMPTPStyrene->AddConstProperty("SCINTILLATIONYIELD",10./keV);
431 
433 }
static const double cm
Definition: G4SIunits.hh:118
static const double MeV
Definition: G4SIunits.hh:211
G4MaterialPropertiesTable * fMPTPStyrene
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:58
void AddConstProperty(const char *key, G4double PropertyValue)
void ReinitializeGeometry(G4bool destroyFirst=false, G4bool prop=true)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
G4MaterialPropertiesTable * fLXe_mt
static const double keV
Definition: G4SIunits.hh:213
static const double mm
Definition: G4SIunits.hh:114
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:446
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetDimensions()

void LXeDetectorConstruction::SetDimensions ( G4ThreeVector  dims)

Definition at line 358 of file LXeDetectorConstruction.cc.

358  {
359  this->fScint_x=dims[0];
360  this->fScint_y=dims[1];
361  this->fScint_z=dims[2];
363 }
void ReinitializeGeometry(G4bool destroyFirst=false, G4bool prop=true)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetHousingReflectivity()

void LXeDetectorConstruction::SetHousingReflectivity ( G4double  r)

Definition at line 444 of file LXeDetectorConstruction.cc.

444  {
445  fRefl=r;
447 }
void ReinitializeGeometry(G4bool destroyFirst=false, G4bool prop=true)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetHousingThickness()

void LXeDetectorConstruction::SetHousingThickness ( G4double  d_mtl)

Definition at line 367 of file LXeDetectorConstruction.cc.

367  {
368  this->fD_mtl=d_mtl;
370 }
void ReinitializeGeometry(G4bool destroyFirst=false, G4bool prop=true)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetMainScintYield()

void LXeDetectorConstruction::SetMainScintYield ( G4double  y)

Definition at line 472 of file LXeDetectorConstruction.cc.

472  {
473  fLXe_mt->AddConstProperty("SCINTILLATIONYIELD",y/MeV);
474 }
static const double MeV
Definition: G4SIunits.hh:211
Double_t y
void AddConstProperty(const char *key, G4double PropertyValue)
G4MaterialPropertiesTable * fLXe_mt
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetMainVolumeOn()

void LXeDetectorConstruction::SetMainVolumeOn ( G4bool  b)

Definition at line 458 of file LXeDetectorConstruction.cc.

458  {
461 }
void ReinitializeGeometry(G4bool destroyFirst=false, G4bool prop=true)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetNFibers()

void LXeDetectorConstruction::SetNFibers ( G4int  n)

Definition at line 465 of file LXeDetectorConstruction.cc.

465  {
466  fNfibers=n;
468 }
Char_t n[5]
void ReinitializeGeometry(G4bool destroyFirst=false, G4bool prop=true)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetNX()

void LXeDetectorConstruction::SetNX ( G4int  nx)

Definition at line 374 of file LXeDetectorConstruction.cc.

374  {
375  this->fNx=nx;
377 }
void ReinitializeGeometry(G4bool destroyFirst=false, G4bool prop=true)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetNY()

void LXeDetectorConstruction::SetNY ( G4int  ny)

Definition at line 381 of file LXeDetectorConstruction.cc.

381  {
382  this->fNy=ny;
384 }
void ReinitializeGeometry(G4bool destroyFirst=false, G4bool prop=true)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetNZ()

void LXeDetectorConstruction::SetNZ ( G4int  nz)

Definition at line 388 of file LXeDetectorConstruction.cc.

388  {
389  this->fNz=nz;
391 }
void ReinitializeGeometry(G4bool destroyFirst=false, G4bool prop=true)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetPMTRadius()

void LXeDetectorConstruction::SetPMTRadius ( G4double  outerRadius_pmt)

Definition at line 395 of file LXeDetectorConstruction.cc.

395  {
396  this->fOuterRadius_pmt=outerRadius_pmt;
398 }
void ReinitializeGeometry(G4bool destroyFirst=false, G4bool prop=true)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetSphereOn()

void LXeDetectorConstruction::SetSphereOn ( G4bool  b)

Definition at line 437 of file LXeDetectorConstruction.cc.

437  {
438  fSphereOn=b;
440 }
void ReinitializeGeometry(G4bool destroyFirst=false, G4bool prop=true)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetWLSScintYield()

void LXeDetectorConstruction::SetWLSScintYield ( G4double  y)

Definition at line 478 of file LXeDetectorConstruction.cc.

478  {
479  fMPTPStyrene->AddConstProperty("SCINTILLATIONYIELD",y/MeV);
480 }
static const double MeV
Definition: G4SIunits.hh:211
G4MaterialPropertiesTable * fMPTPStyrene
Double_t y
void AddConstProperty(const char *key, G4double PropertyValue)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetWLSSlabOn()

void LXeDetectorConstruction::SetWLSSlabOn ( G4bool  b)

Definition at line 451 of file LXeDetectorConstruction.cc.

451  {
452  fWLSslab=b;
454 }
void ReinitializeGeometry(G4bool destroyFirst=false, G4bool prop=true)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fAir

G4Material* LXeDetectorConstruction::fAir
private

Definition at line 117 of file LXeDetectorConstruction.hh.

◆ fAl

G4Material* LXeDetectorConstruction::fAl
private

Definition at line 114 of file LXeDetectorConstruction.hh.

◆ fC

G4Element* LXeDetectorConstruction::fC
private

Definition at line 119 of file LXeDetectorConstruction.hh.

◆ fD_mtl

G4double LXeDetectorConstruction::fD_mtl
private

Definition at line 131 of file LXeDetectorConstruction.hh.

◆ fDetectorMessenger

LXeDetectorMessenger* LXeDetectorConstruction::fDetectorMessenger
private

Definition at line 106 of file LXeDetectorConstruction.hh.

◆ fExperimentalHall_box

G4Box* LXeDetectorConstruction::fExperimentalHall_box
private

Definition at line 108 of file LXeDetectorConstruction.hh.

◆ fExperimentalHall_log

G4LogicalVolume* LXeDetectorConstruction::fExperimentalHall_log
private

Definition at line 109 of file LXeDetectorConstruction.hh.

◆ fExperimentalHall_phys

G4VPhysicalVolume* LXeDetectorConstruction::fExperimentalHall_phys
private

Definition at line 110 of file LXeDetectorConstruction.hh.

◆ fGlass

G4Material* LXeDetectorConstruction::fGlass
private

Definition at line 121 of file LXeDetectorConstruction.hh.

◆ fH

G4Element* LXeDetectorConstruction::fH
private

Definition at line 120 of file LXeDetectorConstruction.hh.

◆ fLXe

G4Material* LXeDetectorConstruction::fLXe
private

Definition at line 113 of file LXeDetectorConstruction.hh.

◆ fLXe_mt

G4MaterialPropertiesTable* LXeDetectorConstruction::fLXe_mt
private

Definition at line 145 of file LXeDetectorConstruction.hh.

◆ fMainVolume

LXeMainVolume* LXeDetectorConstruction::fMainVolume
private

Definition at line 143 of file LXeDetectorConstruction.hh.

◆ fMainVolumeOn

G4bool LXeDetectorConstruction::fMainVolumeOn
private

Definition at line 140 of file LXeDetectorConstruction.hh.

◆ fMPTPStyrene

G4MaterialPropertiesTable* LXeDetectorConstruction::fMPTPStyrene
private

Definition at line 146 of file LXeDetectorConstruction.hh.

◆ fN

G4Element* LXeDetectorConstruction::fN
private

Definition at line 115 of file LXeDetectorConstruction.hh.

◆ fNfibers

G4int LXeDetectorConstruction::fNfibers
private

Definition at line 136 of file LXeDetectorConstruction.hh.

◆ fNx

G4int LXeDetectorConstruction::fNx
private

Definition at line 132 of file LXeDetectorConstruction.hh.

◆ fNy

G4int LXeDetectorConstruction::fNy
private

Definition at line 133 of file LXeDetectorConstruction.hh.

◆ fNz

G4int LXeDetectorConstruction::fNz
private

Definition at line 134 of file LXeDetectorConstruction.hh.

◆ fO

G4Element* LXeDetectorConstruction::fO
private

Definition at line 116 of file LXeDetectorConstruction.hh.

◆ fOuterRadius_pmt

G4double LXeDetectorConstruction::fOuterRadius_pmt
private

Definition at line 135 of file LXeDetectorConstruction.hh.

◆ fPethylene1

G4Material* LXeDetectorConstruction::fPethylene1
private

Definition at line 124 of file LXeDetectorConstruction.hh.

◆ fPethylene2

G4Material* LXeDetectorConstruction::fPethylene2
private

Definition at line 125 of file LXeDetectorConstruction.hh.

◆ fPMMA

G4Material* LXeDetectorConstruction::fPMMA
private

Definition at line 123 of file LXeDetectorConstruction.hh.

◆ fPmt_SD

G4Cache<LXePMTSD*> LXeDetectorConstruction::fPmt_SD
private

Definition at line 150 of file LXeDetectorConstruction.hh.

◆ fPstyrene

G4Material* LXeDetectorConstruction::fPstyrene
private

Definition at line 122 of file LXeDetectorConstruction.hh.

◆ fRefl

G4double LXeDetectorConstruction::fRefl
private

Definition at line 138 of file LXeDetectorConstruction.hh.

◆ fScint_SD

G4Cache<LXeScintSD*> LXeDetectorConstruction::fScint_SD
private

Definition at line 149 of file LXeDetectorConstruction.hh.

◆ fScint_x

G4double LXeDetectorConstruction::fScint_x
private

Definition at line 128 of file LXeDetectorConstruction.hh.

◆ fScint_y

G4double LXeDetectorConstruction::fScint_y
private

Definition at line 129 of file LXeDetectorConstruction.hh.

◆ fScint_z

G4double LXeDetectorConstruction::fScint_z
private

Definition at line 130 of file LXeDetectorConstruction.hh.

◆ fSlab_z

G4double LXeDetectorConstruction::fSlab_z
private

Definition at line 141 of file LXeDetectorConstruction.hh.

◆ fSphereOn

G4bool LXeDetectorConstruction::fSphereOn = true
staticprivate

Definition at line 137 of file LXeDetectorConstruction.hh.

◆ fVacuum

G4Material* LXeDetectorConstruction::fVacuum
private

Definition at line 118 of file LXeDetectorConstruction.hh.

◆ fWLSslab

G4bool LXeDetectorConstruction::fWLSslab
private

Definition at line 139 of file LXeDetectorConstruction.hh.


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