BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSMaterials.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#include "BDSDebug.hh"
20#include "BDSException.hh"
21#include "BDSMaterials.hh"
22#include "BDSParser.hh"
23#include "BDSUtilities.hh"
24#include "BDSWarning.hh"
25
26#include "G4MaterialTable.hh"
27#include "G4String.hh"
28#include "G4NistManager.hh"
29#include "G4Version.hh"
30
31#include <iomanip>
32#include <list>
33#include <map>
34#include <set>
35#include <vector>
36
38
40{
41 if (!instance)
42 {instance = new BDSMaterials();}
43 return instance;
44}
45
46BDSMaterials::BDSMaterials()
47{
48#ifdef BDSDEBUG
49 G4cout << "BDSMaterials: G4 predefined units: " << G4endl;
50 G4cout << "g= " << CLHEP::g << G4endl;
51 G4cout << "m= " << CLHEP::m << G4endl;
52 G4cout << "mole= " << CLHEP::mole << G4endl;
53 G4cout << "kelvin= " << CLHEP::kelvin << G4endl;
54#endif
55
65}
66
68{
69 //default Geant4 temperature = 293.15 K (NTP_Temperature)
70 //default Geant4 pressure = 1atm
71#if G4VERSION_NUMBER < 1011
72 G4double NTP_Temperature = 293.15;
73#endif
74
75 // solid materials
76 // metals
77 // standard single element metals
78 AddExistingMaterialAlias("Al", "aluminium");
79 AddExistingMaterialAlias("Be", "beryllium");
80 AddExistingMaterialAlias("C", "carbon");
81 AddExistingMaterialAlias("Cr", "chromium");
82 AddExistingMaterialAlias("Cu", "copper");
83 AddExistingMaterialAlias("Fe", "iron");
84 AddExistingMaterialAlias("Pb", "lead");
85 AddExistingMaterialAlias("Mg", "magnesium");
86 AddExistingMaterialAlias("Ni", "nickel");
87 AddExistingMaterialAlias("Si", "silicon");
88 AddExistingMaterialAlias("Ti", "titanium");
89 AddExistingMaterialAlias("W", "tungsten");
90 AddExistingMaterialAlias("U", "uranium");
91 AddExistingMaterialAlias("V", "vanadium");
92 AddExistingMaterialAlias("Zn", "zinc");
93
94 // special forms
95 std::list<int> singleElement = {1};
96 AddMaterial("graphite", 2.265, kStateSolid, NTP_Temperature, 1, {"C"}, singleElement);
97 AddMaterial("graphitefoam", 0.61, kStateSolid, NTP_Temperature, 1, {"C"}, singleElement);
98 AddMaterial("solidhydrogen",8.96, kStateSolid, NTP_Temperature, 1, {"H"}, singleElement);
99 AddMaterial("solidnitrogen",8.96, kStateSolid, NTP_Temperature, 1, {"N"}, singleElement);
100 AddMaterial("solidoxygen", 8.96, kStateSolid, NTP_Temperature, 1, {"O"}, singleElement);
101 AddMaterial("weightiron", 7.87, kStateSolid, NTP_Temperature, 1, {"Fe"}, singleElement);
102
103 // composites and alloys
104 AddMaterial("leadtungstate" ,
105 8.27,
106 kStateSolid, NTP_Temperature, 1,
107 {"Pb","W" ,"O" },
108 std::list<int>{1, 1, 4});
109
110 AddMaterial("smco",
111 8.4,
112 kStateSolid, 300, 1,
113 {"Sm","Co"},
114 std::list<double>{0.338, 0.662});
115
116 // Titanium alloy (BDS collimator material SLAC-TN-03-069 p25) deacon 15 Jun 2007
117 AddMaterial("titaniumalloy",
118 4.48,
119 kStateSolid, 300, 1,
120 {"V" ,"Al","Ti"},
121 std::list<double>{0.025, 0.03, 0.945});
122
123 // Carbon Steel (shell of cryomodule). LDeacon 21 Feb 2006
124 AddMaterial("carbonsteel",
125 7.87,
126 kStateSolid, 100, 1,
127 {"C","Mn","P","S","Fe"},
128 std::list<double>{0.0017, 0.0045, 0.0004, 0.0005, 0.9929});
129
130 // Copper Alloy 17410 (C17410) "Beryllium Copper"
131 AddMaterial("berylliumcopper",
132 8.8,
133 kStateSolid, 300, 1,
134 {"Cu", "Be", "Co", "Al", "Fe", "Ni"},
135 std::list<double>{0.991, 0.0031, 0.00500, 0.0004, 0.0003, 0.0002});
136
137 // Stainless Steel 316L
138 AddMaterial("stainlesssteel",
139 8.0,
140 kStateSolid, 295, 1,
141 {"C","Mn","Si","P","S","Cr","Mo","Ni","N","Fe"},
142 std::list<double>{0.0003, 0.02, 0.0075, 0.00045, 0.0003, 0.17, 0.025, 0.12, 0.001, 0.65545});
143
144 // Stainless Steel AISI code 304L (low-carbon) @ 300K
145 AddMaterial("stainless_steel_304L",
146 8.02,
147 kStateSolid, 300, 1,
148 {"Fe", "Cr", "Ni", "Mn", "Si", "P", "S", "C"},
149 std::list<double>{0.67145, 0.185, 0.1125, 0.02, 0.01, 0.00045, 0.0003, 0.0003});
150
151 // Stainless Steel AISI code 304L (low-carbon) @ 87K
152 AddMaterial("stainless_steel_304L_87K",
153 8.02,
154 kStateSolid, 87, 1,
155 {"Fe", "Cr", "Ni", "Mn", "Si", "P", "S", "C"},
156 std::list<double>{0.67145, 0.185, 0.1125, 0.02, 0.01, 0.00045, 0.0003, 0.0003});
157
158 // Stainless Steel AISI code 304L (low-carbon) @ 2K
159 AddMaterial("stainless_steel_304L_2K",
160 8.02,
161 kStateSolid, 2, 1,
162 {"Fe", "Cr", "Ni", "Mn", "Si", "P", "S", "C"},
163 std::list<double>{0.67145, 0.185, 0.1125, 0.02, 0.01, 0.00045, 0.0003, 0.0003});
164
165 // Stainless Steel AISI code 316LN
166 // (Type 316, low carbon, nitrogen-enhanced) @ 300K
167 AddMaterial("stainless_steel_316LN",
168 8.03,
169 kStateSolid, 300, 1,
170 {"Fe", "Cr", "Ni", "Mo", "Mn", "Si", "Ti", "N", "Nb",
171 "Cu", "Co", "P", "C", "S", "Ta", "B"},
172 std::list<double>{0.65093, 0.1700, 0.12000, 0.02500, 0.0200, 0.00750,
173 0.00150, 0.0014, 0.00100, 0.00100, 0.0005, 0.00045,
174 0.00030, 0.0003, 0.00010, 0.00002});
175
176 // Stainless Steel AISI code 316LN
177 // (Type 316, low-carbon nitrogen-enhanced) @ 87K
178 AddMaterial("stainless_steel_316LN_87K",
179 8.03,
180 kStateSolid, 87, 1,
181 {"Fe", "Cr", "Ni", "Mo", "Mn", "Si", "Ti", "N",
182 "Nb", "Cu", "Co", "P", "C", "S", "Ta", "B"},
183 std::list<double>{0.65093, 0.1700, 0.12000, 0.02500, 0.0200,
184 0.00750, 0.00150, 0.0014, 0.00100, 0.00100,
185 0.0005, 0.00045, 0.00030, 0.0003, 0.00010,
186 0.00002});
187
188 // Stainless Steel AISI code 316LN
189 // (Type 316, low-carbon nitrogen-enhanced) @ 2K
190 AddMaterial("stainless_steel_316LN_2K",
191 8.03,
192 kStateSolid, 2, 1,
193 {"Fe", "Cr", "Ni", "Mo", "Mn", "Si", "Ti", "N",
194 "Nb", "Cu", "Co", "P", "C", "S", "Ta", "B"},
195 std::list<double>{0.65093, 0.1700, 0.12000, 0.02500, 0.0200,
196 0.00750, 0.00150, 0.0014, 0.00100, 0.00100,
197 0.0005, 0.00045, 0.00030, 0.0003, 0.00010,
198 0.00002});
199
200
201 // Mild Steel
202 AddMaterial("mild_steel", 8.000, kStateSolid, 295, 1,
203 {"C", "Mn", "Si", "Fe"},
204 std::list<double>{0.002, 0.005, 0.0015, 0.99150});
205
206 // Pure tungsten is not typically used, but instead part of "heavy
207 // alloy." I think "heavy" in the sense that the tungsten makes up almost
208 // all of the composition, and tungsten is a very dense metal.
209 AddMaterial("tungsten_heavy_alloy",
210 18.5,
211 kStateSolid, 87, 1,
212 {"W", "Ni", "Fe"},
213 std::list<double>{0.97, 0.02, 0.01});
214
215 // For HL-LHC Collimation - Inermet170
216 AddMaterial("inermet170",
217 17.0,
218 kStateSolid, 300, 1,
219 {"W", "Ni", "Cu"},
220 std::list<double>{0.90, 0.05, 0.05});
221
222 // For HL-LHC Collimation - Inermet176
223 AddMaterial("inermet176",
224 17.6,
225 kStateSolid, 300, 1,
226 {"W", "Ni", "Cu"},
227 std::list<double>{0.925, 0.0375, 0.0375});
228
229 // For HL-LHC Collimation - Inermet180
230 AddMaterial("inermet180",
231 18.0,
232 kStateSolid, 300, 1,
233 {"W", "Ni", "Cu"},
234 std::list<double>{0.95, 0.025, 0.025});
235}
236
238{
239 // niobium at 2K
240 AddMaterial("niobium_2k", 8.57 , kStateSolid, 2, 1, {"Nb"}, std::list<int>{1});
241 AddExistingMaterialAlias("niobium_2k", "nb_2k"); // alias
242 // niobium titanium at 4K
243 AddMaterial("nbti_4k", 5.6 , kStateSolid, 4, 1, {"Nb","Ti"}, std::list<int>{1,1});
244}
245
247{
248 // boron Nickel (absorber)
249 AddMaterial("bn5000",
250 1.925,
251 kStateSolid, 300, 1,
252 {"B","Ni","O","Ca","Si"},
253 std::list<double>{0.383249242, 0.472071387, 0.0366276887, 0.0228923054, 0.0851593762});
254
255 // calcium carbonate (calcite)
256 AddMaterial("calciumCarbonate",
257 2.711,
258 kStateSolid, 300, 1,
259 {"Ca","O","C"},
260 std::list<int>{1,3,1});
261
262 // clay
263 AddMaterial("clay",
264 1.746,
265 kStateSolid, 300, 1,
266 {"Al","O","Si","H"},
267 std::list<int>{1,9,2,4});
268
269 AddMaterial("concrete",
270 2.3,
271 kStateSolid, 300, 1, {"Si","O","H","Ca","Al","Fe"},
272 std::list<double>{0.227915, 0.60541, 0.09972, 0.04986, 0.014245, 0.00285});
273
274 // we make the material manually here so we have a pointer to it. note density units must
275 // be used here.
276 G4Material* tmpMaterial = new G4Material("fusedsilica",
277 1.032*CLHEP::g/CLHEP::cm3, // density
278 2,
279 kStateSolid);
280 tmpMaterial->AddElement(GetElement("O") , 2);
281 tmpMaterial->AddElement(GetElement("Si"), 1);
282 const G4int FusedSilica_NUMENTRIES = 3; //Number of entries in the material properties table
283 G4double FusedSilica_RIND[FusedSilica_NUMENTRIES]={1.49,1.49,1.49};
284 G4double FusedSilica_AbsLength[FusedSilica_NUMENTRIES]={420.*CLHEP::cm,420.*CLHEP::cm,420.*CLHEP::cm};
285 G4double FusedSilica_Energy[FusedSilica_NUMENTRIES] = {2.0*CLHEP::eV,7.0*CLHEP::eV,7.14*CLHEP::eV};
286 G4MaterialPropertiesTable* fsMaterialPropertiesTable = CreatePropertiesTable();
287 fsMaterialPropertiesTable->AddProperty("ABSLENGTH",FusedSilica_Energy,FusedSilica_AbsLength,FusedSilica_NUMENTRIES);
288 fsMaterialPropertiesTable->AddProperty("RINDEX",FusedSilica_Energy,FusedSilica_RIND,FusedSilica_NUMENTRIES);
289 tmpMaterial->SetMaterialPropertiesTable(fsMaterialPropertiesTable);
290 AddMaterial(tmpMaterial, "fusedsilica");
291
292 // perspex
293 AddMaterial("perspex",
294 1.18,
295 kStateSolid,
296 300,
297 1,
298 {"C","O","H"},
299 std::list<double>{0.59984,0.31961,0.08055});
300
301 // invar - Temperature 2 kelvin. LDeacon 6th Feburary 2006
302 AddMaterial("invar" ,
303 8.1 ,
304 kStateSolid, 2, 1,
305 {"Ni","Fe"},
306 std::list<double>{0.35,0.65});
307
308 // kapton polyimide film
309 AddMaterial("kapton",
310 1.42,
311 kStateSolid, 295, 1,
312 {"H","C","N","O"},
313 std::list<double>{0.026362,0.691133,0.073270,0.209235});
314
315 // n-bk7
316 tmpMaterial = new G4Material("n-bk7", 1.032*CLHEP::g/CLHEP::cm3, 2, kStateSolid);
317 tmpMaterial->AddElement(GetElement("O") , 2);
318 tmpMaterial->AddElement(GetElement("Si"), 1);
319 const G4int N_Bk7_NUMENTRIES = 3; //Number of entries in the material properties table
320 G4double N_Bk7_RIND[N_Bk7_NUMENTRIES]={1.51680,1.51680,1.51680};
321 G4double N_Bk7_AbsLength[N_Bk7_NUMENTRIES]={420.*CLHEP::cm,420.*CLHEP::cm,420.*CLHEP::cm};
322 G4double N_Bk7_Energy[N_Bk7_NUMENTRIES] = {2.0*CLHEP::eV,7.0*CLHEP::eV,7.14*CLHEP::eV};
323 G4MaterialPropertiesTable* nbk7MaterialPropertiesTable = CreatePropertiesTable();
324 nbk7MaterialPropertiesTable->AddProperty("ABSLENGTH",N_Bk7_Energy,N_Bk7_AbsLength,N_Bk7_NUMENTRIES);
325 nbk7MaterialPropertiesTable->AddProperty("RINDEX",N_Bk7_Energy,N_Bk7_RIND,N_Bk7_NUMENTRIES);
326 tmpMaterial->SetMaterialPropertiesTable(nbk7MaterialPropertiesTable);
327 AddMaterial(tmpMaterial, "n-bk7");
328
329 // quartz
330 AddMaterial("quartz", 2.655, kStateSolid, 300, 1, {"Si","O"}, std::list<int>{1,2});
331
332 // marl
333 AddMaterial("marl",
334 1.474,
335 kStateSolid, 300, 1,
336 {"clay","calciumCarbonate"},
337 std::list<double>{0.5,0.5});
338
339 // clayousMarl
340 AddMaterial("clayousMarl",
341 1.555,
342 kStateSolid, 300, 1,
343 {"clay","calciumCarbonate"},
344 std::list<double>{0.65,0.35});
345
346 // limousMarl
347 AddMaterial("limousMarl",
348 1.392,
349 kStateSolid, 300, 1,
350 {"clay", "calciumCarbonate"},
351 std::list<double>{0.35,0.65});
352
353 // "standard" soil (dry)
354 AddMaterial("soil",
355 1.9,
356 kStateSolid, 300, 1,
357 {"Si","O","H","Al"},
358 std::list<double>{0.33377483443708611, 0.57218543046357617, 0.022516556291390728, 0.071523178807947022});
359
360 // epoxy Resin components
361 // the main component of epoxy resin commonly used to insulate magnet coils.
362 AddMaterial("aralditef", 1.175, kStateSolid, 300, 1, {"C","H","O"},std::list<int>{12,18,4});
363
364 // a hardener for the epoxy resin
365 AddMaterial("hy906", 1.225, kStateSolid, 300, 1, {"C","H","O"},std::list<int>{10,5,3});
366
367 // an accelerator for epoxy resin
368 AddMaterial("dy061", 1.025, kStateSolid, 300, 1, {"C","H","O","N"},std::list<int>{15,25,1,3});
369
370 // material type 3 from CERN 81-05, "The Selection and Properties of Epoxide Resins
371 // Used for the Insulation of Magnet Systems in Radiation Environments".
372 AddMaterial("epoxyresin3",
373 1.20,
374 kStateSolid, 300, 1,
375 {"aralditef","hy906","dy061"},
376 std::list<double>{0.497512,0.497512,0.004976});
377
378 // cellulose
379 tmpMaterial = G4NistManager::Instance()->FindOrBuildMaterial("G4_CELLULOSE_CELLOPHANE",true,true);
380 const G4int Cellulose_NUMENTRIES = 3; //Number of entries in the material properties table
381 G4double Cellulose_RIND[Cellulose_NUMENTRIES] = {1.532,1.532,1.532};//Assume constant refractive index.
382 G4double Cellulose_Energy[Cellulose_NUMENTRIES] = {2.0*CLHEP::eV,7.0*CLHEP::eV,7.14*CLHEP::eV}; //The energies.
383 G4MaterialPropertiesTable* celluloseMaterialPropertiesTable = CreatePropertiesTable();
384 celluloseMaterialPropertiesTable->AddProperty("RINDEX",Cellulose_Energy, Cellulose_RIND, Cellulose_NUMENTRIES);
385 tmpMaterial->SetMaterialPropertiesTable(celluloseMaterialPropertiesTable);
386 AddMaterial(tmpMaterial, "cellulose");
387
388 // polyurethane
389 AddMaterial("polyurethane",
390 1.05,
391 kStateSolid, NTP_Temperature, 1,
392 {"C","H","N","O"},
393 std::list<int>{6,10,2,4});
394
395 // RCH 1000 - Ultra high molecular weight polyethylene [PE-UHMW]
396 // at 4K for LHC dipoles
397 AddMaterial("rch1000_4k",
398 0.925,
399 kStateSolid, 4, 1,
400 {"C", "H"},
401 std::list<int>{2,4});
402}
403
405{
406 //YAG
407 G4Material* tmpMaterial = new G4Material("yag", 4.56*CLHEP::g/CLHEP::cm3, 3);
408 tmpMaterial->AddElement(GetElement("Y") , 3);
409 tmpMaterial->AddElement(GetElement("Al"), 5);
410 tmpMaterial->AddElement(GetElement("O") , 12);
411 G4double birks = 0.08*CLHEP::mm/CLHEP::MeV;
412 tmpMaterial->GetIonisation()->SetBirksConstant(birks);
413 G4MaterialPropertiesTable* mpt_YAG = CreatePropertiesTable();
414 // const int nEntries3=60;
415 const G4int nEntries = 9;
416 G4double PhotonEnergyYAG[nEntries];
417 G4double dNEntries2=(G4double)nEntries;
418 G4double energyMin=1.91*CLHEP::eV;
419 G4double energyMax=2.76*CLHEP::eV;
420 G4double deltaEnergy=(energyMax-energyMin)/(dNEntries2-1.0);
421 G4double energy=energyMin;
422 for (G4int i=0; i<nEntries; energy += deltaEnergy, i++)
423 {PhotonEnergyYAG[i]=energy;}
424
425 G4double RefractiveIndexYAG[nEntries] = //Approximately correct, but check for different wavelengths
426 { 1.82, 1.82, 1.82, 1.82, 1.82, 1.82, 1.82,
427 1.82, 1.82 };
428
429 mpt_YAG->AddProperty("RINDEX",PhotonEnergyYAG, RefractiveIndexYAG, nEntries);
430#if G4VERSION_NUMBER < 1079
431 G4double scintFastYAG[nEntries] = //Approximately correct
432 { 0, 0.25, 2.0, 14.0, 13.0, 7.0, 4.0, 2.0, 0.0 };
433 // All of these parameters are deprecated in V11 onwards and will cause a crash
434 mpt_YAG->AddProperty("FASTCOMPONENT",PhotonEnergyYAG, scintFastYAG, nEntries)->SetSpline(true);
435 mpt_YAG->AddConstProperty("FASTTIMECONSTANT",70.*CLHEP::ns); //Approximately correct
436 mpt_YAG->AddConstProperty("YIELDRATIO",1.0);
437#endif
438 mpt_YAG->AddConstProperty("SCINTILLATIONYIELD",8000./CLHEP::MeV); //Approximately correct
439 mpt_YAG->AddConstProperty("RESOLUTIONSCALE",2.0); //Check this
440 tmpMaterial->SetMaterialPropertiesTable(mpt_YAG);
441 AddMaterial(tmpMaterial, "yag");
442
443 //UPS-923A - see http://www.amcrys-h.com/
444 //Define the material properties (copy from NIST table of materials).
445 G4NistManager* nistManager = G4NistManager::Instance();
446 G4Material* polystyrene = nistManager->FindOrBuildMaterial("G4_POLYSTYRENE",true,true);
447 tmpMaterial = new G4Material("ups923a",polystyrene->GetDensity(),1);
448 tmpMaterial->AddMaterial(polystyrene,1);
449 tmpMaterial->SetName("ups923a");
450 std::vector<G4double> ups923a_PhotonEnergy = {
451 3.35, 3.31, 3.28, 3.26, 3.25, 3.23, 3.23,
452 3.22, 3.21, 3.19, 3.18, 3.17, 3.16, 3.15,
453 3.14, 3.14, 3.13, 3.11, 3.1, 3.09, 3.09,
454 3.08, 3.07, 3.04, 3.02, 3.02, 3.01, 2.99,
455 2.98, 2.97, 2.97, 2.95, 2.95, 2.93, 2.93,
456 2.92, 2.92, 2.91, 2.89, 2.88, 2.87, 2.86,
457 2.85, 2.83, 2.81, 2.8, 2.79, 2.78, 2.76,
458 2.74, 2.72, 2.71, 2.68, 2.66, 2.64, 2.62,
459 2.61, 2.58, 2.55, 2.53, 2.5, 2.48, 2.46,
460 2.44, 2.41, 2.38, 2.35 };
461 std::reverse(ups923a_PhotonEnergy.begin(), ups923a_PhotonEnergy.end());
462
463 // AUG 21 - these were previously just one number as a const property but they should be non-const
464 // which requires vs energy numbers - so just use arrays the same shape
465 std::vector<G4double> ups923a_RINDEX(ups923a_PhotonEnergy.size(), 1.52);
466 std::vector<G4double> ups923a_ABSLENGTH(ups923a_PhotonEnergy.size(), 1*CLHEP::m);
467
468 G4MaterialPropertiesTable* ups923a_mt = CreatePropertiesTable();
469#if G4VERSION_NUMBER < 1070
470 ups923a_mt->AddProperty("RINDEX", ups923a_PhotonEnergy.data(), ups923a_RINDEX.data(), (int)ups923a_PhotonEnergy.size());
471 ups923a_mt->AddProperty("ABSLENGTH", ups923a_PhotonEnergy.data(), ups923a_ABSLENGTH.data(), (int)ups923a_PhotonEnergy.size());
472#else
473 ups923a_mt->AddProperty("RINDEX", ups923a_PhotonEnergy, ups923a_RINDEX);
474 ups923a_mt->AddProperty("ABSLENGTH", ups923a_PhotonEnergy, ups923a_ABSLENGTH);
475#endif
476 //Birk's constant
477 birks = (0.014/1.06)*CLHEP::cm/CLHEP::MeV;
478 tmpMaterial->GetIonisation()->SetBirksConstant(birks);
479#if G4VERSION_NUMBER < 1079
480 const G4int ups923a_numentries = 67;
481 G4double ups923a_emission[ups923a_numentries] = {
482 0, 0.04, 0.11, 0.2, 0.3, 0.4, 0.52,
483 0.62, 0.67, 0.68, 0.67, 0.62, 0.53, 0.48,
484 0.44, 0.42, 0.4, 0.41, 0.42, 0.51, 0.46,
485 0.57, 0.67, 0.78, 0.91, 0.93, 0.95, 0.96,
486 0.94, 0.91, 0.85, 0.76, 0.67, 0.61, 0.57,
487 0.55, 0.52, 0.51, 0.52, 0.54, 0.57, 0.58,
488 0.6, 0.6, 0.59, 0.58, 0.55, 0.48, 0.42,
489 0.37, 0.33, 0.31, 0.29, 0.28, 0.26, 0.24,
490 0.2, 0.17, 0.12, 0.09, 0.08, 0.07,
491 0.06, 0.04, 0.02, 0.01, 0.01 };
492 ups923a_mt->AddConstProperty("FASTTIMECONSTANT",3.3*CLHEP::ns);
493 ups923a_mt->AddProperty("FASTCOMPONENT",ups923a_PhotonEnergy.data(), ups923a_emission, ups923a_numentries)->SetSpline(true);
494 ups923a_mt->AddConstProperty("YIELDRATIO",1.0);
495#endif
496 ups923a_mt->AddConstProperty("RESOLUTIONSCALE",2.0); //Check this
497 G4double scintYieldAnthracene=14200; //Anthracene yield per 1 CLHEP::MeV
498 G4double scintYieldUPS923A=scintYieldAnthracene*0.60;//60% of anthracene
499 ups923a_mt->AddConstProperty("SCINTILLATIONYIELD",scintYieldUPS923A/CLHEP::MeV);
500
501 tmpMaterial->SetMaterialPropertiesTable(ups923a_mt);
502 AddMaterial(tmpMaterial, "ups923a");
503
504 // PET
505 G4double pet_density=1.4*CLHEP::g/CLHEP::cm3;
506 G4int pet_nelements=3;
507 G4State pet_state=kStateSolid;
508 tmpMaterial= new G4Material("pet",
509 pet_density,
510 pet_nelements,
511 pet_state);
512 tmpMaterial->AddElement(nistManager->FindOrBuildElement("C",true),10);
513 tmpMaterial->AddElement(nistManager->FindOrBuildElement("H",true),8);
514 tmpMaterial->AddElement(nistManager->FindOrBuildElement("O",true),4);
515 const G4int Pet_NUMENTRIES = 3; //Number of entries in the material properties table
516 G4double Pet_RIND[Pet_NUMENTRIES] = {1.570,1.570,1.570};//Assume constant refractive index.
517 G4double Pet_Energy[Pet_NUMENTRIES] = {2.0*CLHEP::eV,7.0*CLHEP::eV,7.14*CLHEP::eV}; //The energies.
518 G4MaterialPropertiesTable* petMaterialPropertiesTable = CreatePropertiesTable();
519 petMaterialPropertiesTable->AddProperty("RINDEX",Pet_Energy, Pet_RIND, Pet_NUMENTRIES);
520 tmpMaterial->SetMaterialPropertiesTable(petMaterialPropertiesTable);
521 AddMaterial(tmpMaterial, "pet");
522
523 // Opaque PET (Dacron)
524 tmpMaterial= new G4Material("pet_opaque",
525 pet_density,
526 pet_nelements,
527 pet_state);
528 tmpMaterial->AddElement(nistManager->FindOrBuildElement("C",true),10);
529 tmpMaterial->AddElement(nistManager->FindOrBuildElement("H",true),8);
530 tmpMaterial->AddElement(nistManager->FindOrBuildElement("O",true),4);
531 const G4int Pet_Opaque_NUMENTRIES = 3; //Number of entries in the material properties table
532 G4double Pet_Opaque_RIND[Pet_Opaque_NUMENTRIES] = {1.570,1.570,1.570};//Assume constant refractive index.
533 G4double Pet_Opaque_Energy[Pet_Opaque_NUMENTRIES] = {2.0*CLHEP::eV,7.0*CLHEP::eV,7.14*CLHEP::eV}; //The energies.
534 G4double Pet_Opaque_abslen[]={1*CLHEP::um, 1*CLHEP::um, 1*CLHEP::um};
535 G4MaterialPropertiesTable* pet_opaqueMaterialPropertiesTable = CreatePropertiesTable();
536 pet_opaqueMaterialPropertiesTable->AddProperty("RINDEX",Pet_Opaque_Energy, Pet_Opaque_RIND, Pet_Opaque_NUMENTRIES);
537 pet_opaqueMaterialPropertiesTable->AddProperty("ABSLENGTH",Pet_Opaque_Energy, Pet_Opaque_abslen, Pet_Opaque_NUMENTRIES);
538 tmpMaterial->SetMaterialPropertiesTable(pet_opaqueMaterialPropertiesTable);
539 AddMaterial(tmpMaterial, "pet_opaque");
540
541 // Gadolinium oxysulphate Gd_2 O_2 S
542 G4Material* GOS = nistManager->FindOrBuildMaterial("G4_GADOLINIUM_OXYSULFIDE",true,true);
543
544 // Ganolinium oxysulphate in a polyurethane elastomer (lanex)
545 G4double fill_factor=0.5;
546 G4double lanex_density=fill_factor*GOS->GetDensity()+(1-fill_factor)*GetMaterial("polyurethane")->GetDensity();
547 G4double gos_fraction_by_mass=fill_factor*GOS->GetDensity()/lanex_density;
548 G4double pur_fraction_by_mass=1-gos_fraction_by_mass;
549
550 tmpMaterial = new G4Material("lanex", lanex_density, 2);
551 tmpMaterial->AddMaterial(GOS, gos_fraction_by_mass);
552 tmpMaterial->AddMaterial(GetMaterial("polyurethane"), pur_fraction_by_mass);
553
554 G4MaterialPropertiesTable* mptLanex = CreatePropertiesTable();
555 const G4int nentLanex=2;
556 G4double rindex=1.50;//(1.82+1.50)/2.0;
557 G4double energytab[]={2.239*CLHEP::eV, 2.241*CLHEP::eV};
558 G4double rindextab[]={rindex, rindex};
559 G4double abslen[]={7*CLHEP::mm, 7*CLHEP::mm};
560 mptLanex->AddProperty("RINDEX",energytab, rindextab, nentLanex); //Average refractive index of bulk material
561 mptLanex->AddProperty("ABSLENGTH", energytab, abslen, nentLanex);
562#if G4VERSION_NUMBER < 1079
563 mptLanex->AddConstProperty("MIEHG", 60.3e-3*CLHEP::mm); // interface changed in V11
564 G4double emitspec[]={1.0, 1.0};
565 mptLanex->AddProperty("FASTCOMPONENT",energytab, emitspec, nentLanex);
566 mptLanex->AddConstProperty("FASTTIMECONSTANT", 1.*CLHEP::ns);
567#endif
568 G4double scintScalingFactor=1;
569 mptLanex->AddConstProperty("SCINTILLATIONYIELD",7.8e4/CLHEP::MeV);
570 mptLanex->AddConstProperty("RESOLUTIONSCALE",1.0);
571 mptLanex->AddConstProperty("MIEHG_FORWARD", 0.91);
572 mptLanex->AddConstProperty("MIEHG_BACKWARD", 0.91);
573 mptLanex->AddConstProperty("MIEHG_FORWARD_RATIO", 1.0);
574 tmpMaterial->SetMaterialPropertiesTable(mptLanex);
575 AddMaterial(tmpMaterial, "lanex");
576
577 //Ganolinium oxysulphate in a polyurethane elastomer (lanex) - version 2
578 tmpMaterial = new G4Material("lanex2", lanex_density, 2);
579 tmpMaterial->AddMaterial(GOS, gos_fraction_by_mass);
580 tmpMaterial->AddMaterial(GetMaterial("polyurethane"), pur_fraction_by_mass);
581 G4MaterialPropertiesTable* mptLanex2 = CreatePropertiesTable();
582 mptLanex2->AddProperty("RINDEX",energytab, rindextab, nentLanex); //Average refractive index of bulk material
583 mptLanex2->AddProperty("ABSLENGTH", energytab, abslen, nentLanex);
584#if G4VERSION_NUMBER < 1079
585 mptLanex2->AddConstProperty("MIEHG", 60.3e-3*CLHEP::mm);
586 mptLanex2->AddProperty("FASTCOMPONENT",energytab, emitspec, nentLanex);
587 mptLanex2->AddConstProperty("FASTTIMECONSTANT", 1.*CLHEP::ns);
588#endif
589 mptLanex2->AddConstProperty("SCINTILLATIONYIELD",8.9e4/CLHEP::MeV);
590 mptLanex2->AddConstProperty("RESOLUTIONSCALE",1.0);
591 mptLanex2->AddConstProperty("MIEHG_FORWARD", 0.91);
592 mptLanex2->AddConstProperty("MIEHG_BACKWARD", 0.91);
593 mptLanex2->AddConstProperty("MIEHG_FORWARD_RATIO", 0.5);
594 tmpMaterial->SetMaterialPropertiesTable(mptLanex);
595 AddMaterial(tmpMaterial, "lanex2");
596
597 //gos_lanex - GOS with the bulk optical transport properties of lanex particles suspended in an elastomer but the atomic, density and scintillation properties of GOS
598 G4double gos_lanex_density=GOS->GetDensity();
599 tmpMaterial = new G4Material("gos_lanex", gos_lanex_density, 1);
600 tmpMaterial->AddMaterial(GOS, 1.0);
601 G4MaterialPropertiesTable* mptGOSLanex = CreatePropertiesTable();
602 const G4int nentGOSLanex=2;
603 G4double rindexGOSLanex=1.50;
604 G4double energyGOSLanexTab[]={2.239*CLHEP::eV, 2.241*CLHEP::eV};
605 G4double rindexGOSLanexTab[]={rindexGOSLanex, rindexGOSLanex};
606 G4double abslenGOSLanex[]={7*CLHEP::mm, 7*CLHEP::mm};
607 G4double gosLanexMiehgForward=0.911;
608 G4double gosLanexMiehgBackward=0.911;
609 G4double gosLanexMiehgForwardRatio=0.5;
610#if G4VERSION_NUMBER < 1079
611 G4double mieHgTimeConst=1.0*CLHEP::ns;
612 G4double emitspecGOSLanex[]={1.0, 1.0};
613 G4double mieScatteringLengthGOSLanex=60.3*CLHEP::um;
614 mptGOSLanex->AddProperty("FASTCOMPONENT",energyGOSLanexTab, emitspecGOSLanex, nentGOSLanex);
615 mptGOSLanex->AddConstProperty("FASTTIMECONSTANT", mieHgTimeConst);
616 mptGOSLanex->AddConstProperty("YIELDRATIO", 1.0);
617 mptGOSLanex->AddConstProperty("MIEHG", mieScatteringLengthGOSLanex);
618#endif
619 mptGOSLanex->AddConstProperty("SCINTILLATIONYIELD",8.9e4/CLHEP::MeV); //Intrinisic scintilation yield of GOS
620 mptGOSLanex->AddConstProperty("RESOLUTIONSCALE", 1.0);
621 mptGOSLanex->AddConstProperty("MIEHG_FORWARD", gosLanexMiehgForward);
622 mptGOSLanex->AddConstProperty("MIEHG_BACKWARD", gosLanexMiehgBackward);
623 mptGOSLanex->AddConstProperty("MIEHG_FORWARD_RATIO", gosLanexMiehgForwardRatio);
624 mptGOSLanex->AddProperty("RINDEX",energyGOSLanexTab, rindexGOSLanexTab, nentGOSLanex); //Average refractive index of bulk material
625 mptGOSLanex->AddProperty("ABSLENGTH", energyGOSLanexTab, abslenGOSLanex, nentGOSLanex);
626 tmpMaterial->SetMaterialPropertiesTable(mptGOSLanex);
627 AddMaterial(tmpMaterial, "gos_lanex");
628
629 //Same as gos_lanex but refractive index = 1
630 tmpMaterial = new G4Material("gos_ri1", gos_lanex_density, 1);
631 tmpMaterial->AddMaterial(GOS, 1.0);
632 G4MaterialPropertiesTable* mptGOSLanexRi1 = CreatePropertiesTable();
633 G4double rindexGOSLanexRi1Tab[]={1.0, 1.0};
634#if G4VERSION_NUMBER < 1079
635 mptGOSLanexRi1->AddProperty("FASTCOMPONENT",energyGOSLanexTab, emitspecGOSLanex, nentGOSLanex);
636 mptGOSLanexRi1->AddConstProperty("FASTTIMECONSTANT", mieHgTimeConst);
637 mptGOSLanexRi1->AddConstProperty("YIELDRATIO", 1.0);
638 mptGOSLanexRi1->AddConstProperty("MIEHG", mieScatteringLengthGOSLanex);
639#endif
640 mptGOSLanexRi1->AddConstProperty("SCINTILLATIONYIELD",8.9e4/CLHEP::MeV); //Intrinisic scintilation yield of GOS
641 mptGOSLanexRi1->AddConstProperty("RESOLUTIONSCALE", 1.0);
642 mptGOSLanexRi1->AddConstProperty("MIEHG_FORWARD", gosLanexMiehgForward);
643 mptGOSLanexRi1->AddConstProperty("MIEHG_BACKWARD", gosLanexMiehgBackward);
644 mptGOSLanexRi1->AddConstProperty("MIEHG_FORWARD_RATIO", gosLanexMiehgForwardRatio);
645 mptGOSLanexRi1->AddProperty("RINDEX",energyGOSLanexTab, rindexGOSLanexRi1Tab, nentGOSLanex); //Average refractive index of bulk material
646 mptGOSLanexRi1->AddProperty("ABSLENGTH", energyGOSLanexTab, abslenGOSLanex, nentGOSLanex);
647 tmpMaterial->SetMaterialPropertiesTable(mptGOSLanexRi1);
648 AddMaterial(tmpMaterial, "gos_ri1");
649
650 //pet_lanex - PET with the bulk optical transport properties of lanex particles suspended in an elastomer but the atomic, density and scintillation properties of PET
651 G4double pet_lanex_density=GetMaterial("polyurethane")->GetDensity();
652 tmpMaterial = new G4Material("pet_lanex", pet_lanex_density, 1);
653 tmpMaterial->AddMaterial(GetMaterial("polyurethane"), 1.0);
654 G4MaterialPropertiesTable* mptPETLanex = CreatePropertiesTable();
655 mptPETLanex->AddConstProperty("MIEHG_FORWARD", gosLanexMiehgForward);
656 mptPETLanex->AddConstProperty("MIEHG_BACKWARD", gosLanexMiehgBackward);
657 mptPETLanex->AddConstProperty("MIEHG_FORWARD_RATIO", gosLanexMiehgForwardRatio);
658#if G4VERSION_NUMBER < 1079
659 mptPETLanex->AddConstProperty("MIEHG", mieScatteringLengthGOSLanex);
660#endif
661 mptPETLanex->AddProperty("RINDEX",energyGOSLanexTab, rindexGOSLanexTab, nentGOSLanex); //Average refractive index of bulk material
662 mptPETLanex->AddProperty("ABSLENGTH", energyGOSLanexTab, abslenGOSLanex, nentGOSLanex);
663 tmpMaterial->SetMaterialPropertiesTable(mptPETLanex);
664 AddMaterial(tmpMaterial, "pet_lanex");
665
666 //Medex (larger grained lanex)
667 // G4double medex_fill_factor=0.5;
668 G4double medex_density=fill_factor*GOS->GetDensity()+(1-fill_factor)*GetMaterial("polyurethane")->GetDensity();
669 G4double medex_gos_fraction_by_mass=fill_factor*GOS->GetDensity()/medex_density;
670 G4double medex_pur_fraction_by_mass=1-medex_gos_fraction_by_mass;
671 tmpMaterial = new G4Material("medex", medex_density, 2);
672 tmpMaterial->AddMaterial(GOS, medex_gos_fraction_by_mass);
673 tmpMaterial->AddMaterial(GetMaterial("polyurethane"), medex_pur_fraction_by_mass);
674 G4MaterialPropertiesTable* mptMedex = CreatePropertiesTable();
675 const G4int nentMedex=2;
676 // G4double medexRindex=(1.82+1.50)/2.0;
677 // G4double medexEnergytab[]={2.239*CLHEP::eV, 2.241*CLHEP::eV};
678 G4double medexRindextab[]={rindex, rindex};
679
680 G4double medexAbslen[]={7*CLHEP::mm, 7*CLHEP::mm};
681 mptMedex->AddProperty("RINDEX",energytab, medexRindextab, nentMedex); //Average refractive index of bulk material
682 mptMedex->AddProperty("ABSLENGTH", energytab, medexAbslen, nentMedex);
683#if G4VERSION_NUMBER < 1079
684 mptMedex->AddConstProperty("MIEHG", 230e-3*CLHEP::mm);
685 G4double medexEmitspec[]={1.0, 1.0};
686 mptMedex->AddProperty("FASTCOMPONENT",energytab, medexEmitspec, nentMedex);
687 mptMedex->AddConstProperty("FASTTIMECONSTANT", 1.*CLHEP::ns);
688#endif
689 mptMedex->AddConstProperty("SCINTILLATIONYIELD",scintScalingFactor*2.94e4/CLHEP::MeV);
690 mptMedex->AddConstProperty("RESOLUTIONSCALE",1.0);
691 mptMedex->AddConstProperty("MIEHG_FORWARD", 0.93);
692 mptMedex->AddConstProperty("MIEHG_BACKWARD", 0.93);
693 mptMedex->AddConstProperty("MIEHG_FORWARD_RATIO", 1.0);
694 tmpMaterial->SetMaterialPropertiesTable(mptMedex);
695 AddMaterial(tmpMaterial, "medex");
696
697 // carbon fiber
698 G4double frac_graph=0.5;
699 G4double frac_poly=0.5;
700 G4Material* graph=GetMaterial("G4_GRAPHITE");
701 G4Material* poly=GetMaterial("G4_POLYACRYLONITRILE");
702 G4double dens_graph=graph->GetDensity();
703 G4double dens_poly=poly->GetDensity();
704
705 G4double dens_cf =frac_graph*dens_graph+frac_poly*dens_poly;
706 G4double frac_graph_bw=frac_graph*dens_graph/dens_cf;
707 G4double frac_poly_bw=frac_poly*dens_poly/dens_cf;
708 tmpMaterial = new G4Material("carbonfiber", dens_cf, 2);
709 tmpMaterial->AddMaterial(graph, frac_graph_bw);
710 tmpMaterial->AddMaterial(poly, frac_poly_bw);
711 G4MaterialPropertiesTable* mptCarbonfiber = CreatePropertiesTable();
712 const G4int nentCarbonfiber=2;
713 G4double energytab_cf[]={2.239*CLHEP::eV, 2.241*CLHEP::eV};
714 G4double carbonfiberRindextab[]={2.6, 2.6};
715 G4double carbonfiberAbslen[]={2.1*CLHEP::um, 2.1*CLHEP::um};
716 mptCarbonfiber->AddProperty("RINDEX",energytab_cf, carbonfiberRindextab, nentCarbonfiber); //Average refractive index of bulk material
717 mptCarbonfiber->AddProperty("ABSLENGTH", energytab_cf, carbonfiberAbslen, nentCarbonfiber);
718 tmpMaterial->SetMaterialPropertiesTable(mptCarbonfiber);
719 AddMaterial(tmpMaterial, "carbonfiber");
720}
721
723{
724 AddMaterial("lhcconcrete",
725 2.42,
726 kStateSolid, 300, 1,
727 {"H","C","O","Na","Mg","Al","Si","K","Ca","Fe","P","S",
728 "Ti","Mn","Zn","Zr","Ba","Pb","Sr","Eu"},
729 std::list<double>{0.59785345499811 *CLHEP::perCent,
730 5.59989402848226 *CLHEP::perCent,
731 49.1111702720319 *CLHEP::perCent,
732 0.45137935852357 *CLHEP::perCent,
733 0.66062806777291 *CLHEP::perCent,
734 2.05561946276849 *CLHEP::perCent,
735 18.7995018924154 *CLHEP::perCent,
736 0.65365311079793 *CLHEP::perCent,
737 20.0191229406116 *CLHEP::perCent,
738 1.11400027114647 *CLHEP::perCent,
739 0.04782827639985 *CLHEP::perCent,
740 0.01195706909996 *CLHEP::perCent,
741 0.3457585814739 *CLHEP::perCent,
742 0.03856154784738 *CLHEP::perCent,
743 0.02401378044242 *CLHEP::perCent,
744 0.00737352594498 *CLHEP::perCent,
745 0.01783596140744 *CLHEP::perCent,
746 0.04623400051985 *CLHEP::perCent,
747 0.39757254757374 *CLHEP::perCent,
748 4.184974185E-05 *CLHEP::perCent});
749
750 // LHC Rock
751 AddMaterial("lhc_rock",
752 2,
753 kStateSolid, 300, 1,
754 {"H", "Na", "Si", "Fe", "C", "Mg", "K", "O", "Al", "Ca"},
755 std::list<double>{0.006, 0.01, 0.2, 0.014, 0.03, 0.005, 0.01, 0.5, 0.03, 0.195});
756
757 // Superconducting components of LHC magnet elements
758 // Definitions taken from FLUKA
759
760 std::list<int> singleElement = {1};
761 // Liquid helium at 1.9K
762 AddMaterial("lhe_1.9k",
763 0.1472,
764 kStateLiquid, 1.9, 1,
765 {"He"}, singleElement);
766
767 // Niobium @ 87K
768 AddMaterial("nb_87k",
769 8.902 ,
770 kStateSolid , 87 , 1,
771 {"Nb"}, singleElement);
772
773 // Titanium @ 87K
774 AddMaterial("ti_87k",
775 4.54,
776 kStateSolid, 87 , 1,
777 {"Ti"}, singleElement);
778
779 // superconductor NbTi with Ti = 47% by weight
780 AddMaterial("nbti_87k",
781 6.0471,
782 kStateSolid , 87 , 1,
783 {"nb_87k","ti_87k"},
784 std::list<double>{0.53,0.47});
785
786 // copper at 4 Kelvin
787 AddMaterial("cu_4k",
788 8.96,
789 kStateSolid, 4, 1,
790 {"Cu"}, singleElement);
791
792 // copper at 2 Kelvin
793 AddMaterial("cu_2k",
794 8.96,
795 kStateSolid, 2, 1,
796 {"Cu"}, singleElement);
797
798 // naked superconductor NbTi wire with Cu/SC volume ratio (>= 4.0 and <4.8)
799 AddMaterial("nbti.1",
800 8.4206,
801 kStateSolid, 4, 1,
802 {"nbti_87k","cu_4k"},
803 std::list<double>{1.0/5.4, 4.4/5.4});
804
805 // For HL-LHC Collimation - Copper Diamond - TBC - component fractions
806 AddMaterial("copperdiamond",
807 5.3,
808 kStateSolid, 300, 1,
809 {"Cu", "C"},
810 std::list<double>{0.32, 0.68});
811 AddExistingMaterialAlias("copperdiamond", "cucd");
812
813 // For HL-LHC Collimation - Copper Diamond
814 AddMaterial("molybdenumcarbide",
815 8.9,
816 kStateSolid, 300, 1,
817 {"Mo", "C"},
818 std::list<int>{1, 1});
819 AddExistingMaterialAlias("molybdenumcarbide", "mogr");
820}
821
823{
824 AddMaterial("liquidhelium", 0.12498, kStateLiquid, 4.15, 1, {"He"}, std::list<int>{1});
825
826 // water
827 AddExistingMaterialAlias("G4_WATER", "water"); // use NIST water
828 G4Material* water = GetMaterial("water"); // get the material pointer back
829
830 // add refractive index to G4_WATER material for cherenkov radiation to work
831 const G4int nEntries = 9;
832 G4MaterialPropertiesTable* waterProperties = CreatePropertiesTable();
833 G4double photonEnergy[nEntries];
834 G4double dNEntries = (G4double)nEntries;
835 G4double energyMin = 1.*CLHEP::eV;
836 G4double energyMax = 3.*CLHEP::eV;
837 G4double deltaEnergy = (energyMax-energyMin)/(dNEntries-1.0);
838 G4double energy = energyMin;
839 for (G4int i = 0; i < nEntries; energy += deltaEnergy, i++)
840 {photonEnergy[i] = energy;}
841 G4double refractiveIndex[nEntries] = {1.325, 1.325, 1.326, 1.327, 1.328, 1.33, 1.333, 1.336, 1.343};
842 waterProperties->AddProperty("RINDEX", photonEnergy, refractiveIndex, nEntries);
843 water->SetMaterialPropertiesTable(waterProperties);
844}
845
847{
848 AddExistingMaterialAlias("N", "nitrogen");
849
850 // air
851 AddExistingMaterialAlias("G4_AIR", "air");
852 G4Material* g4AirMaterial = GetMaterial("air");
853 const G4int airNEntries = 3;
854 // source: NPL Tables of Physical & Chemical Constants. Refractive indices at
855 // different energies.
856 G4double airRefractiveIndex[airNEntries] = {1.000292,1.000292,1.000292};
857 G4double airEnergy[airNEntries] = {2.0*CLHEP::eV,7.0*CLHEP::eV,7.14*CLHEP::eV};
858 G4MaterialPropertiesTable* airProperties = CreatePropertiesTable();
859 airProperties->AddProperty("RINDEX", airEnergy, airRefractiveIndex, airNEntries);
860 g4AirMaterial->SetMaterialPropertiesTable(airProperties);
861
862 // old way - manually constructed air
863 // room temperature temperature and standard pressure
864 G4double temperature = 300*CLHEP::kelvin;
865 G4double pressure = 1.0*CLHEP::atmosphere;
866 G4double airDensity = 0.001225; // g/cm3
867 G4Material* tmpMaterial = new G4Material("airbdsim",
868 airDensity*CLHEP::g/CLHEP::cm3,
869 2,
870 kStateGas,
871 temperature,
872 pressure);
873 tmpMaterial->AddElement(GetElement("O"), 0.2);
874 tmpMaterial->AddElement(GetElement("N"), 0.8);
875 tmpMaterial->SetMaterialPropertiesTable(airProperties);
876 AddMaterial(tmpMaterial, "airbdsim");
877
878 // Carbon monoxide
879 G4double coDensity = 0.001145; // g/cm3
880 AddMaterial("carbonmonoxide",
881 coDensity,
882 kStateGas,
883 temperature,
884 pressure,
885 {"C","O"},
886 std::list<int>{1,1});
887
888 // Carbon monoxide beam pipe gas
889 G4double bp_pressure = 0.0133e-9 * (CLHEP::bar/CLHEP::atmosphere); //10 nTorr pressure
890 AddMaterial("bp_carbonmonoxide",
891 coDensity,
892 kStateGas,
893 temperature,
894 bp_pressure,
895 {"C","O"},
896 std::list<int>{1,1});
897}
898
900{
901 //Awake plasma - rubidium at density of 7e14 atoms/cm3
902 // G4double numberDensity = 7.0e14/CLHEP::cm3;
903 G4double a = 85.4678*CLHEP::g/CLHEP::mole;
904 G4double density = 1e-7 * CLHEP::g/CLHEP::cm3;
905 G4Material* tmpMaterial = new G4Material("awakeplasma", 37., a, density);
906 AddMaterial(tmpMaterial, "awakeplasma");
907}
908
910{
911 // Default vacuum (same composition as residual vacuum in warm sections of LHC).
912 G4double vacpressure;
913 // check if parser is initialised (in case list materials is asked by command line)
915 {vacpressure=BDSParser::Instance()->GetOptions().vacuumPressure*CLHEP::bar;}
916 else
917 {vacpressure=1e-12*CLHEP::bar;}
918 G4double temperature = 300*CLHEP::kelvin;
919 G4double density = (CLHEP::STP_Temperature/temperature) * (vacpressure/(1.*CLHEP::atmosphere)) * 29*CLHEP::g/(22.4*1.e-3*CLHEP::m3);
920
921 G4Material* tmpMaterial = new G4Material("vacuum",
922 density,
923 3, kStateGas,
924 temperature, vacpressure);
925 tmpMaterial->AddElement(GetElement("H"), 0.482);
926 tmpMaterial->AddElement(GetElement("C"), 0.221);
927 tmpMaterial->AddElement(GetElement("O"), 0.297);
928 AddMaterial(tmpMaterial, "vacuum");
929
930 const G4int Vac_NUMENTRIES = 3;
931 // Assume refractive index = 1 in a vacuum.
932 G4double Vac_RIND[Vac_NUMENTRIES] = {1.000,1.000,1.000};
933 G4double Vac_Energy[Vac_NUMENTRIES] = {2.0*CLHEP::eV,7.0*CLHEP::eV,7.14*CLHEP::eV};
934 G4MaterialPropertiesTable* vacMaterialPropertiesTable = CreatePropertiesTable();
935 vacMaterialPropertiesTable->AddProperty("RINDEX",Vac_Energy, Vac_RIND, Vac_NUMENTRIES);
936 tmpMaterial->SetMaterialPropertiesTable(vacMaterialPropertiesTable);
937
938 // copy regular vacuum as named 'laservac' to detect where laserwire is.
939 // G4Material doesn't have a copy constructor but it has a sort of copy
940 // constructor.
941 G4Material* regularVacuum = GetMaterial("vacuum");
942 G4Material* laservac = new G4Material("laservac",
943 regularVacuum->GetDensity(),
944 regularVacuum,
945 kStateGas,
946 regularVacuum->GetTemperature(),
947 regularVacuum->GetPressure());
948 AddMaterial(laservac, "laservac");
949}
950
951void BDSMaterials::AddMaterial(G4Material* material, G4String name)
952{
953 name = BDS::LowerCase(name);
954 if (materialNames.count(name) == 0) // not defined already
955 {
956 materials[name] = material;
957 materialNames.insert(name);
958 }
959 else
960 {throw BDSException(__METHOD_NAME__, "Material \"" + name + "\" already exists");}
961}
962
963void BDSMaterials::AddExistingMaterialAlias(const G4String &existingMaterialName,
964 G4String alias)
965{
966 alias = BDS::LowerCase(alias);
967 G4Material* material = GetMaterial(existingMaterialName);
968 aliases[alias] = material; // store in lower case as that's how we search
969 aliasNames.insert(alias);
970
971 // cache aliases for a given material for feedback
972 auto search = materialToAliases.find(material);
973 if (search == materialToAliases.end())
974 {materialToAliases[material] = {};}
975 materialToAliases[material].push_back(alias);
976}
977
978void BDSMaterials::AddMaterial(G4String name,
979 G4double Z,
980 G4double A,
981 G4double density,
982 G4State state,
983 G4double temperature,
984 G4double pressure)
985{
986 // convention: material name in small letters (to be able to find materials regardless of capitalisation)
987 name = BDS::LowerCase(name);
988 DensityCheck(density, name);
989
990 G4Material* tmpMaterial = new G4Material(name,
991 Z,
992 A*CLHEP::g/CLHEP::mole,
993 density*CLHEP::g/CLHEP::cm3,
994 state,
995 temperature*CLHEP::kelvin,
996 pressure*CLHEP::atmosphere);
997 AddMaterial(tmpMaterial, name);
998}
999
1000template <typename Type>
1001void BDSMaterials::AddMaterial(G4String name,
1002 G4double density,
1003 G4State state,
1004 G4double temperature,
1005 G4double pressure,
1006 const std::list<G4String>& components,
1007 const std::list<Type>& componentFractions)
1008{
1009 name = BDS::LowerCase(name);
1010 DensityCheck(density, name);
1011
1012 G4Material* tmpMaterial = new G4Material(name,
1013 density*CLHEP::g/CLHEP::cm3,
1014 (G4int)components.size(),
1015 state,
1016 temperature*CLHEP::kelvin,
1017 pressure*CLHEP::atmosphere);
1018 std::list<G4String>::const_iterator sIter;
1019 typename std::list<Type>::const_iterator dIter;
1020 for (sIter = components.begin(), dIter = componentFractions.begin();
1021 sIter != components.end();
1022 ++sIter, ++dIter)
1023 {
1024#ifdef BDSDEBUG
1025 G4cout << "BDSMaterials::AddMaterial: " << *sIter << G4endl;
1026#endif
1027 G4Element* element = CheckElement(*sIter);
1028 if (element)
1029 {tmpMaterial->AddElement(element, (*dIter));}
1030 else
1031 {tmpMaterial->AddMaterial(GetMaterial(*sIter), (*dIter));}
1032 }
1033 AddMaterial(tmpMaterial, name);
1034}
1035
1036G4Material* BDSMaterials::GetMaterial(G4String material) const
1037{
1038 if (material.empty())
1039 {throw BDSException(__METHOD_NAME__, "empty material name");}
1040 G4String materialOriginal = material;
1041 // for short names we assume they're elements so we prefix with G4_ and
1042 // get them from NIST
1043 G4String nistString ("G4_");
1044 if (material.length() <= 2)
1045#if G4VERSION_NUMBER > 1099
1046 {material = nistString + material;}
1047#else
1048 {material.prepend(nistString);}
1049#endif
1050
1051 G4String start = material.substr(0,3);
1052 if (nistString == start)
1053 {
1054#ifdef BDSDEBUG
1055 G4cout << "Using NIST material " << material << G4endl;
1056#endif
1057 G4Material* mat = G4NistManager::Instance()->FindOrBuildMaterial(material, true, true);
1058 if (!mat)
1059 {throw BDSException(__METHOD_NAME__, "\"" + material + "\" could not be found by NIST.");}
1060 return mat;
1061 }
1062 else
1063 {
1064 // find material regardless of capitalisation
1065 material = BDS::LowerCase(material);
1066
1067 // search predefined + parser materials
1068 if (materialNames.count(material) == 1)
1069 {return materials.at(material);}
1070 else if (aliasNames.count(material) == 1)
1071 {return aliases.at(material);}
1072 else
1073 {
1074 // search aliases
1075 if (externalMaterialNames.count(material) == 1)
1076 {return externalMaterials.at(material);}
1077 else
1078 {// can't find it -> warn and exit
1079 ListMaterials();
1080 throw BDSException(__METHOD_NAME__, "\"" + materialOriginal + "\" is unknown.");
1081 }
1082 }
1083 }
1084}
1085
1086void BDSMaterials::CheckForConflictingMaterialsAfterLoad(const G4String& geometryFileName,
1087 const G4String& componentName) const
1088{
1089 G4MaterialTable* table = G4Material::GetMaterialTable();
1090
1091 G4bool check = false;
1092 std::map<G4String, G4int> nameCount;
1093 for (const auto& mat : *table)
1094 {
1095 G4String name = mat->GetName();
1096 auto search = nameCount.find(name);
1097 if (search == nameCount.end())
1098 {nameCount[name] = 0;}
1099 nameCount[name] += 1;
1100 check = check || nameCount[name] > 1; // flag for efficiency - don't loop again if we don't need to
1101 }
1102
1103 if (check)
1104 {
1105 for (const auto& nc: nameCount)
1106 {
1107 if (nc.second > 1)
1108 {
1109 G4String msg = "the material \"" + nc.first + "\" has been defined more\n";
1110 msg += "than once now and will cause the wrong material to be used for any future usages.\n";
1111 msg += "Error caused by GDML file \""+geometryFileName+"\" in component \""+componentName+"\"\n";
1112 msg += "This is most likely due to a material name in the GDML file conflicting with a predefined one in BDSIM.";
1113 if (BDS::StartsWith(BDS::LowerCase(nc.first), "g4_"))
1114 {BDS::Warning(msg);} // only warn for a NIST one as ultimately it won't be that different even if redefined from reloading
1115 else
1116 {throw BDSException(__METHOD_NAME__, msg);}
1117 }
1118 }
1119 }
1120}
1121
1122
1123void BDSMaterials::CacheMaterialsFromGDML(const std::map<G4String, G4Material*>& materialsGDML)
1124{
1125 // this function is really just concerning GDML material caching and not
1126 // all external format name caching...
1127 for (const auto& kv : materialsGDML)
1128 {
1129 G4String nameLower = BDS::LowerCase(kv.first);
1130
1131 // We have to allow names that are in bdsim materials as the user could have
1132 // used just a bdsim material by name. However, we can't know if that was one
1133 // that was constructed uniquely with the same name. For this, we need to check
1134 // the material table in the above function.
1135 if (BDS::StartsWith(nameLower, "g4_") || materials.find(nameLower) != materials.end())
1136 {continue;} // a Geant4 material or a BDSIM one
1137
1138 if (externalMaterialNames.count(nameLower) == 1)
1139 {
1140 G4String msg = "the material \""+kv.first+"\" is already defined and has\n";
1141 msg += "already been loaded from a previous GDML file(s) and is ambiguous.\n";
1142 msg += "Please prepend with the BDSIM element used to load the file to be explicit.";
1143 throw BDSException(__METHOD_NAME__, msg);
1144 }
1145
1146 // it's ok for a GDML material to have a material name the same as a BDSIM
1147 // alias as only the real name will conflict
1148
1149 externalMaterials[nameLower] = kv.second;
1150 externalMaterialNames.insert(nameLower);
1151 }
1152}
1153
1154void BDSMaterials::AddElement(G4Element* element, const G4String& symbol)
1155{
1156 if (CheckElement(symbol) != nullptr)
1157 {throw BDSException(__METHOD_NAME__, "Element \"" + symbol + "\" already exists.");}
1158
1159 elements.insert(make_pair(symbol, element));
1160#ifdef BDSDEBUG
1161 G4cout << "New element : " << symbol << G4endl;
1162#endif
1163}
1164
1165void BDSMaterials::AddElement(const G4String& name,
1166 const G4String& symbol,
1167 G4double Z,
1168 G4double A)
1169{
1170 G4Element* tmpElement = new G4Element(name, symbol, Z, A*CLHEP::g/CLHEP::mole);
1171 AddElement(tmpElement, symbol);
1172}
1173
1174void BDSMaterials::DensityCheck(G4double density,
1175 const G4String& materialName) const
1176{
1177 if (density > 1e2)
1178 {// so greater than 100g / cm3, the densest natural material is around 23g/cm3
1179 G4String msg = "material \"" + materialName + "\"has a density higher than 100g/cm3! Perhaps check this!\n";
1180 msg += "Density: " + std::to_string(density) + " g/cm3... proceeding";
1181 BDS::Warning(__METHOD_NAME__, msg);
1182 }
1183}
1184
1185G4Element* BDSMaterials::CheckElement(const G4String& symbol) const
1186{
1187 // first look in defined element list
1188 auto iter = elements.find(symbol);
1189 if(iter != elements.end())
1190 {return (*iter).second;}
1191 else
1192 {
1193 G4Element* element = G4NistManager::Instance()->FindOrBuildElement(symbol, true);
1194 return element;
1195 }
1196}
1197
1198G4Element* BDSMaterials::GetElement(const G4String& symbol) const
1199{
1200 G4Element* element = CheckElement(symbol);
1201 if (!element)
1202 {throw BDSException(__METHOD_NAME__, "Element \"" + symbol + "\" could not be found.");}
1203 return element;
1204}
1205
1206void BDSMaterials::PrintBasicMaterialMassFraction(G4Material* material) const
1207{
1208 // simpler version of geant4 code
1209 for (G4int i = 0; i < (G4int)material->GetNumberOfElements(); i++)
1210 {
1211 G4cout << (*(material->GetElementVector()))[i]->GetName() << "\t "
1212 << (material->GetFractionVector()[i])/CLHEP::perCent << " %" << G4endl;
1213 }
1214}
1215
1217{
1218 auto flagsCache(G4cout.flags());
1219 G4cout << __METHOD_NAME__ << G4endl;
1220 // always print out vacuum composition
1221 G4Material* vacuum = GetMaterial("vacuum");
1222 G4cout<< "\"vacuum\" composition: " << G4endl;
1224 G4cout<< "pressure = " << vacuum->GetPressure()/CLHEP::bar << " bar" << G4endl;
1225 G4cout<< "temperature = " << vacuum->GetTemperature()/CLHEP::kelvin << " K" << G4endl;
1226 G4cout<< "density = " << vacuum->GetDensity()/(CLHEP::g/CLHEP::cm3)<< " g/cm^3" << G4endl << G4endl;
1227
1228 G4cout << "All elements are available with their 1 or 2 letter chemical symbol. ie C or G4_C" << G4endl << G4endl;
1229
1230 if (!elements.empty())
1231 {
1232 G4cout << "Extra defined elements are:" << G4endl;
1233 for (const auto& element : elements)
1234 {G4cout << std::left << std::setw(12) << element.second->GetName() << " - " << element.second->GetSymbol() << G4endl;}
1235 G4cout << G4endl;
1236 }
1237
1238 G4cout << "Defined materials are:" << G4endl;
1239 for (const auto& material : materials)
1240 {
1241 G4cout << material.first;
1242 G4String realName = material.second->GetName();
1243 if (realName != material.first)
1244 {G4cout << " (" << material.second->GetName() << ")" << G4endl;}
1245 else
1246 {G4cout << G4endl;}
1247
1248 auto aliasSearch = materialToAliases.find(material.second);
1249 if (aliasSearch != materialToAliases.end())
1250 {
1251 const auto v = aliasSearch->second;
1252 G4cout << "Aliases: ";
1253 for (const auto& n : v)
1254 {G4cout << "\"" << n << "\" ";}
1255 G4cout << G4endl;
1256 }
1257 }
1258 G4cout << G4endl;
1259 G4cout << "All aliases: (alias, real name)" << G4endl;
1260 for (const auto& aliasMaterial : aliases)
1261 {
1262 G4cout << aliasMaterial.first << " : " << aliasMaterial.second->GetName() << G4endl;}
1263 G4cout << G4endl;
1264 G4cout << "Available NIST materials are:" << G4endl;
1265 G4NistManager::Instance()->ListMaterials("all");
1266 G4cout.flags(flagsCache);
1267}
1268
1269BDSMaterials::~BDSMaterials()
1270{
1271 for (auto& material : materials)
1272 {delete material.second;}
1273 materials.clear();
1274 materialNames.clear();
1275
1276 // we don't delete externalMaterials as any loaded from GDML will be
1277 // deleted by Geant4
1278 externalMaterialNames.clear();
1279
1280 for (auto& element : elements)
1281 {delete element.second;}
1282 elements.clear();
1283
1284 for (G4MaterialPropertiesTable* table : propertiesTables)
1285 {delete table;}
1286
1287 instance = nullptr;
1288}
1289
1291{
1292 // This function uses the list from the parser and prepares
1293 // the necessary materials for this run.
1294#ifdef BDSDEBUG
1295 G4bool debug = true;
1296#else
1297 G4bool debug = false;
1298#endif
1299
1300 // convert the parsed atom list to list of Geant4 G4Elements
1301 if (verbose || debug)
1302 {G4cout << __METHOD_NAME__ << "parsing the atom list..." << G4endl;}
1303 for (const auto& it : BDSParser::Instance()->GetAtoms())
1304 {AddElement(it.name,it.symbol,it.Z,it.A);}
1305 if (verbose || debug)
1306 {G4cout << "size of atom list: "<< BDSParser::Instance()->GetAtoms().size() << G4endl;}
1307
1308 // convert the parsed material list to list of Geant4 G4Materials
1309 if (verbose || debug)
1310 {G4cout << "parsing the material list..."<< G4endl;}
1311 for (auto it : BDSParser::Instance()->GetMaterials())
1312 {
1313 G4State itsState;
1314 if (it.state=="solid")
1315 {itsState = kStateSolid;}
1316 else if (it.state=="liquid")
1317 {itsState = kStateLiquid;}
1318 else if (it.state=="gas")
1319 {itsState = kStateGas;}
1320 else
1321 {
1322 G4cout << "Unknown material state "<< it.state
1323 << ", setting it to default (solid)"
1324 << G4endl;
1325 it.state="solid";
1326 itsState = kStateSolid;
1327 }
1328
1329#ifdef BDSDEBUG
1330 G4cout << "---->adding Material, ";
1331 it.print();
1332#endif
1333
1334 if(it.Z != 0)
1335 {
1336 AddMaterial(it.name,
1337 it.Z,
1338 it.A,
1339 it.density,
1340 itsState,
1341 it.temper,
1342 it.pressure);
1343 }
1344 else if(!(it.components.empty()))
1345 {
1346 std::list<G4String> tempComponents;
1347 for (const auto& jt : it.components)
1348 {tempComponents.emplace_back(G4String(jt));}
1349
1350 if(it.componentsWeights.size()==it.components.size())
1351 {
1352 AddMaterial(it.name,
1353 it.density,
1354 itsState,
1355 it.temper,
1356 it.pressure,
1357 tempComponents,
1358 it.componentsWeights);
1359 }
1360 else if(it.componentsFractions.size()==it.components.size())
1361 {
1362 AddMaterial(it.name,
1363 it.density,
1364 itsState,
1365 it.temper,
1366 it.pressure,
1367 tempComponents,
1368 it.componentsFractions);
1369 }
1370 else
1371 {throw BDSException(__METHOD_NAME__, "Badly defined material - number of components is not equal to number of weights or mass fractions!");}
1372 }
1373 else
1374 {throw BDSException(__METHOD_NAME__, "Badly defined material - need more information!");}
1375 }
1376 if (verbose || debug)
1377 {G4cout << "size of material list: "<< BDSParser::Instance()->GetMaterials().size() << G4endl;}
1378}
1379
1380G4MaterialPropertiesTable* BDSMaterials::CreatePropertiesTable()
1381{
1382 G4MaterialPropertiesTable* table = new G4MaterialPropertiesTable();
1383 propertiesTables.push_back(table);
1384 return table;
1385}
General exception with possible name of object and message.
Definition: BDSException.hh:35
A class for all material definitions known to BDSIM. Additional materials can be added in the parser ...
Definition: BDSMaterials.hh:39
void DefineMetals()
Methods called by constructor.
Definition: BDSMaterials.cc:67
std::set< G4String > aliasNames
Map of other names to existing materials. To avoid double deletion. Also in lower case.
void DefineLHCComponents()
Methods called by constructor.
std::map< G4Material *, std::vector< G4String > > materialToAliases
Keep a vector of possible aliases for every material for print out and user feedback.
void DefineGases()
Methods called by constructor.
void DefineNonMetalSolids()
Methods called by constructor.
void DensityCheck(G4double density, const G4String &materialName) const
Print warning if density suspiciously high. Should be in g/cm3 in G4 units already.
void DefineScintillators()
Methods called by constructor.
void DefineVacuums()
Methods called by constructor.
static BDSMaterials * Instance()
Singleton pattern access.
Definition: BDSMaterials.cc:39
G4MaterialPropertiesTable * CreatePropertiesTable()
Create new properties table and store in vector.
static BDSMaterials * instance
Singleton instance.
std::set< G4String > materialNames
Map of materials, convention name lowercase.
std::map< G4String, G4Element * > elements
Map of elements, no lowercase convention.
void DefineLiquids()
Methods called by constructor.
G4Material * GetMaterial(G4String material) const
Get material by name.
void PrepareRequiredMaterials(G4bool verbose=false)
converts parser material list
void CacheMaterialsFromGDML(const std::map< G4String, G4Material * > &materialsGDML)
Introduce materials loaded from GDML into this instance.
void ListMaterials() const
output available materials
void DefinePlasmas()
Methods called by constructor.
std::set< G4String > externalMaterialNames
Keep hold of externally constructed materials to indicate conflicts to the user but not to delete.
void CheckForConflictingMaterialsAfterLoad(const G4String &geometryFileName, const G4String &componentName) const
void AddElement(G4Element *element, const G4String &symbol)
Add a G4Element.
void AddMaterial(G4Material *material, G4String name)
Add G4Material.
void DefineSuperconductors()
Methods called by constructor.
G4Element * CheckElement(const G4String &symbol) const
Return element if defined (nullptr if not)
void AddExistingMaterialAlias(const G4String &existingMaterialName, G4String alias)
Add alias to a material.
G4Element * GetElement(const G4String &symbol) const
Get element by name.
void PrintBasicMaterialMassFraction(G4Material *material) const
Print mass fractions of consituents of a given material.
std::vector< G4MaterialPropertiesTable * > propertiesTables
Material tables for storing pointers.
static BDSParser * Instance()
Access method.
Definition: BDSParser.cc:28
std::vector< GMAD::Atom > GetAtoms() const
Return the parser list of that object.
Definition: BDSParser.hh:96
static bool IsInitialised()
Returns if parser is initialised.
Definition: BDSParser.cc:49
const GMAD::Options & GetOptions() const
Return options.
Definition: BDSParser.hh:50
std::vector< GMAD::Material > GetMaterials() const
Return the parser list of that object.
Definition: BDSParser.hh:102
G4String LowerCase(const G4String &str)
Utility function to simplify lots of syntax changes for pedantic g4 changes.
G4bool StartsWith(const std::string &expression, const std::string &prefix)
Return true if a string "expression" starts with "prefix".