Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
WLSMaterials.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: WLSMaterials.cc 69561 2013-05-08 12:25:56Z gcosmo $
27 //
28 /// \file optical/wls/src/WLSMaterials.cc
29 /// \brief Implementation of the WLSMaterials class
30 //
31 //
32 #include "WLSMaterials.hh"
33 
34 #include "G4SystemOfUnits.hh"
35 
36 WLSMaterials* WLSMaterials::fInstance = 0;
37 
38 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
39 
40 WLSMaterials::WLSMaterials()
41 {
42  fNistMan = G4NistManager::Instance();
43 
44  fNistMan->SetVerbose(2);
45 
46  CreateMaterials();
47 }
48 
49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50 
52 {
53  delete fPMMA;
54  delete fPethylene;
55  delete fFPethylene;
56  delete fPolystyrene;
57  delete fSilicone;
58 }
59 
60 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61 
63 {
64  if (fInstance == 0)
65  {
66  fInstance = new WLSMaterials();
67  }
68  return fInstance;
69 }
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72 
74 {
75  G4Material* mat = fNistMan->FindOrBuildMaterial(material);
76 
77  if (!mat) mat = G4Material::GetMaterial(material);
78  if (!mat) {
79  std::ostringstream o;
80  o << "Material " << material << " not found!";
81  G4Exception("WLSMaterials::GetMaterial","",
82  FatalException,o.str().c_str());
83  }
84 
85  return mat;
86 }
87 
88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
89 
90 void WLSMaterials::CreateMaterials()
91 {
93  G4int ncomponents;
94  G4double fractionmass;
95  std::vector<G4int> natoms;
96  std::vector<G4double> fractionMass;
97  std::vector<G4String> elements;
98 
99  // Materials Definitions
100  // =====================
101 
102  //--------------------------------------------------
103  // Vacuum
104  //--------------------------------------------------
105 
106  fNistMan->FindOrBuildMaterial("G4_Galactic");
107 
108  //--------------------------------------------------
109  // Air
110  //--------------------------------------------------
111 
112  fAir = fNistMan->FindOrBuildMaterial("G4_AIR");
113 
114  //--------------------------------------------------
115  // WLSfiber PMMA
116  //--------------------------------------------------
117 
118  elements.push_back("C"); natoms.push_back(5);
119  elements.push_back("H"); natoms.push_back(8);
120  elements.push_back("O"); natoms.push_back(2);
121 
122  density = 1.190*g/cm3;
123 
124  fPMMA = fNistMan->
125  ConstructNewMaterial("PMMA", elements, natoms, density);
126 
127  elements.clear();
128  natoms.clear();
129 
130  //--------------------------------------------------
131  // Cladding (polyethylene)
132  //--------------------------------------------------
133 
134  elements.push_back("C"); natoms.push_back(2);
135  elements.push_back("H"); natoms.push_back(4);
136 
137  density = 1.200*g/cm3;
138 
139  fPethylene = fNistMan->
140  ConstructNewMaterial("Pethylene", elements, natoms, density);
141 
142  elements.clear();
143  natoms.clear();
144 
145  //--------------------------------------------------
146  // Double Cladding (fluorinated polyethylene)
147  //--------------------------------------------------
148 
149  elements.push_back("C"); natoms.push_back(2);
150  elements.push_back("H"); natoms.push_back(4);
151 
152  density = 1.400*g/cm3;
153 
154  fFPethylene = fNistMan->
155  ConstructNewMaterial("FPethylene", elements, natoms, density);
156 
157  elements.clear();
158  natoms.clear();
159 
160  //--------------------------------------------------
161  // Polystyrene
162  //--------------------------------------------------
163 
164  elements.push_back("C"); natoms.push_back(8);
165  elements.push_back("H"); natoms.push_back(8);
166 
167  density = 1.050*g/cm3;
168 
169  fPolystyrene = fNistMan->
170  ConstructNewMaterial("Polystyrene", elements, natoms, density);
171 
172  elements.clear();
173  natoms.clear();
174 
175  //--------------------------------------------------
176  // Silicone (Template for Optical Grease)
177  //--------------------------------------------------
178 
179  elements.push_back("C"); natoms.push_back(2);
180  elements.push_back("H"); natoms.push_back(6);
181 
182  density = 1.060*g/cm3;
183 
184  fSilicone = fNistMan->
185  ConstructNewMaterial("Silicone", elements, natoms, density);
186 
187  elements.clear();
188  natoms.clear();
189 
190  //--------------------------------------------------
191  // Aluminium
192  //--------------------------------------------------
193 
194  fNistMan->FindOrBuildMaterial("G4_Al");
195 
196  //--------------------------------------------------
197  // TiO2
198  //--------------------------------------------------
199 
200  elements.push_back("Ti"); natoms.push_back(1);
201  elements.push_back("O"); natoms.push_back(2);
202 
203  density = 4.26*g/cm3;
204 
205  G4Material* TiO2 = fNistMan->
206  ConstructNewMaterial("TiO2", elements, natoms, density);
207 
208  elements.clear();
209  natoms.clear();
210 
211  //--------------------------------------------------
212  // Scintillator Coating - 15% TiO2 and 85% polystyrene by weight.
213  //--------------------------------------------------
214 
215  density = 1.52*g/cm3;
216 
217  fCoating =
218  new G4Material("Coating", density, ncomponents=2);
219 
220  fCoating->AddMaterial(TiO2, fractionmass = 15*perCent);
221  fCoating->AddMaterial(fPolystyrene, fractionmass = 85*perCent);
222 
223  //
224  // ------------ Generate & Add Material Properties Table ------------
225  //
226 
227  const G4int nEntries = 50;
228 
229  G4double photonEnergy[nEntries] =
230  {2.00*eV,2.03*eV,2.06*eV,2.09*eV,2.12*eV,
231  2.15*eV,2.18*eV,2.21*eV,2.24*eV,2.27*eV,
232  2.30*eV,2.33*eV,2.36*eV,2.39*eV,2.42*eV,
233  2.45*eV,2.48*eV,2.51*eV,2.54*eV,2.57*eV,
234  2.60*eV,2.63*eV,2.66*eV,2.69*eV,2.72*eV,
235  2.75*eV,2.78*eV,2.81*eV,2.84*eV,2.87*eV,
236  2.90*eV,2.93*eV,2.96*eV,2.99*eV,3.02*eV,
237  3.05*eV,3.08*eV,3.11*eV,3.14*eV,3.17*eV,
238  3.20*eV,3.23*eV,3.26*eV,3.29*eV,3.32*eV,
239  3.35*eV,3.38*eV,3.41*eV,3.44*eV,3.47*eV};
240 
241  //--------------------------------------------------
242  // Air
243  //--------------------------------------------------
244 
245  G4double refractiveIndex[nEntries] =
246  { 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
247  1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
248  1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
249  1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
250  1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00};
251 
253  mpt->AddProperty("RINDEX", photonEnergy, refractiveIndex, nEntries);
254 
255  fAir->SetMaterialPropertiesTable(mpt);
256 
257  //--------------------------------------------------
258  // PMMA for WLSfibers
259  //--------------------------------------------------
260 
261  G4double refractiveIndexWLSfiber[nEntries] =
262  { 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60,
263  1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60,
264  1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60,
265  1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60,
266  1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60, 1.60};
267 
268  G4double absWLSfiber[nEntries] =
269  {5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,
270  5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,
271  5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,5.40*m,1.10*m,
272  1.10*m,1.10*m,1.10*m,1.10*m,1.10*m,1.10*m, 1.*mm, 1.*mm, 1.*mm, 1.*mm,
273  1.*mm, 1.*mm, 1.*mm, 1.*mm, 1.*mm, 1.*mm, 1.*mm, 1.*mm, 1.*mm, 1.*mm};
274 
275  G4double emissionFib[nEntries] =
276  {0.05, 0.10, 0.30, 0.50, 0.75, 1.00, 1.50, 1.85, 2.30, 2.75,
277  3.25, 3.80, 4.50, 5.20, 6.00, 7.00, 8.50, 9.50, 11.1, 12.4,
278  12.9, 13.0, 12.8, 12.3, 11.1, 11.0, 12.0, 11.0, 17.0, 16.9,
279  15.0, 9.00, 2.50, 1.00, 0.05, 0.00, 0.00, 0.00, 0.00, 0.00,
280  0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00};
281 
282  // Add entries into properties table
284  mptWLSfiber->
285  AddProperty("RINDEX",photonEnergy,refractiveIndexWLSfiber,nEntries);
286  // mptWLSfiber->AddProperty("ABSLENGTH",photonEnergy,absWLSfiber,nEntries);
287  mptWLSfiber->AddProperty("WLSABSLENGTH",photonEnergy,absWLSfiber,nEntries);
288  mptWLSfiber->AddProperty("WLSCOMPONENT",photonEnergy,emissionFib,nEntries);
289  mptWLSfiber->AddConstProperty("WLSTIMECONSTANT", 0.5*ns);
290 
291  fPMMA->SetMaterialPropertiesTable(mptWLSfiber);
292 
293  //--------------------------------------------------
294  // Polyethylene
295  //--------------------------------------------------
296 
297  G4double refractiveIndexClad1[nEntries] =
298  { 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49,
299  1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49,
300  1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49,
301  1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49,
302  1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49, 1.49};
303 
304  G4double absClad[nEntries] =
305  {20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,
306  20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,
307  20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,
308  20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,
309  20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m,20.0*m};
310 
311  // Add entries into properties table
313  mptClad1->AddProperty("RINDEX",photonEnergy,refractiveIndexClad1,nEntries);
314  mptClad1->AddProperty("ABSLENGTH",photonEnergy,absClad,nEntries);
315 
316  fPethylene->SetMaterialPropertiesTable(mptClad1);
317 
318  //--------------------------------------------------
319  // Fluorinated Polyethylene
320  //--------------------------------------------------
321 
322  G4double refractiveIndexClad2[nEntries] =
323  { 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42,
324  1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42,
325  1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42,
326  1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42,
327  1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42, 1.42};
328 
329  // Add entries into properties table
331  mptClad2->AddProperty("RINDEX",photonEnergy,refractiveIndexClad2,nEntries);
332  mptClad2->AddProperty("ABSLENGTH",photonEnergy,absClad,nEntries);
333 
334  fFPethylene->SetMaterialPropertiesTable(mptClad2);
335 
336  //--------------------------------------------------
337  // Silicone
338  //--------------------------------------------------
339 
340  G4double refractiveIndexSilicone[nEntries] =
341  { 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46,
342  1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46,
343  1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46,
344  1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46,
345  1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46, 1.46};
346 
347  // Add entries into properties table
349  mptSilicone->
350  AddProperty("RINDEX",photonEnergy,refractiveIndexSilicone,nEntries);
351  mptSilicone->AddProperty("ABSLENGTH",photonEnergy,absClad,nEntries);
352 
353  fSilicone->SetMaterialPropertiesTable(mptSilicone);
354 
355  //--------------------------------------------------
356  // Polystyrene
357  //--------------------------------------------------
358 
359  G4double refractiveIndexPS[nEntries] =
360  { 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50,
361  1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50,
362  1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50,
363  1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50,
364  1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50, 1.50};
365 
366  G4double absPS[nEntries] =
367  {2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,
368  2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,
369  2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,
370  2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,
371  2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm,2.*cm};
372 
373  G4double scintilFast[nEntries] =
374  {0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
375  0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
376  0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
377  1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0,
378  1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0};
379 
380  // Add entries into properties table
382  mptPolystyrene->AddProperty("RINDEX",photonEnergy,refractiveIndexPS,nEntries);
383  mptPolystyrene->AddProperty("ABSLENGTH",photonEnergy,absPS,nEntries);
384  mptPolystyrene->
385  AddProperty("FASTCOMPONENT",photonEnergy, scintilFast,nEntries);
386  mptPolystyrene->AddConstProperty("SCINTILLATIONYIELD",10./keV);
387  mptPolystyrene->AddConstProperty("RESOLUTIONSCALE",1.0);
388  mptPolystyrene->AddConstProperty("FASTTIMECONSTANT", 10.*ns);
389 
390  fPolystyrene->SetMaterialPropertiesTable(mptPolystyrene);
391 
392  // Set the Birks Constant for the Polystyrene scintillator
393 
394  fPolystyrene->GetIonisation()->SetBirksConstant(0.126*mm/MeV);
395 
396 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
G4Material * GetMaterial(const G4String)
Definition: WLSMaterials.cc:73
void AddMaterial(G4Material *material, G4double fraction)
Definition: G4Material.cc:450
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:578
void SetMaterialPropertiesTable(G4MaterialPropertiesTable *anMPT)
Definition: G4Material.hh:247
void SetBirksConstant(G4double value)
Definition of the WLSMaterials class.
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
G4MaterialPropertyVector * AddProperty(const char *key, G4double *PhotonEnergies, G4double *PropertyValues, G4int NumEntries)
string material
Definition: eplot.py:19
G4double density
Definition: TRTMaterials.hh:39
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
void SetVerbose(G4int)
void AddConstProperty(const char *key, G4double PropertyValue)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual ~WLSMaterials()
Definition: WLSMaterials.cc:51
float perCent
Definition: hepunit.py:239
double G4double
Definition: G4Types.hh:76
static WLSMaterials * GetInstance()
Definition: WLSMaterials.cc:62
#define ns
Definition: xmlparse.cc:597