Geant4-11
Public Member Functions | Private Member Functions | Private Attributes
G4UCNMaterialPropertiesTable Class Reference

#include <G4UCNMaterialPropertiesTable.hh>

Inheritance diagram for G4UCNMaterialPropertiesTable:
G4MaterialPropertiesTable

Public Member Functions

void AddConstProperty (const char *key, G4double propertyValue, G4bool createNewKey=false)
 
void AddConstProperty (const G4String &key, G4double propertyValue, G4bool createNewKey=false)
 
void AddEntry (const char *key, G4double aPhotonEnergy, G4double aPropertyValue)
 
void AddEntry (const G4String &key, G4double aPhotonEnergy, G4double aPropertyValue)
 
G4MaterialPropertyVectorAddProperty (const char *key, G4double *photonEnergies, G4double *propertyValues, G4int numEntries, G4bool createNewKey=false, G4bool spline=false)
 
void AddProperty (const char *key, G4MaterialPropertyVector *opv, G4bool createNewKey=false)
 
void AddProperty (const G4String &key, const G4String &mat)
 
G4MaterialPropertyVectorAddProperty (const G4String &key, const std::vector< G4double > &photonEnergies, const std::vector< G4double > &propertyValues, G4bool createNewKey=false, G4bool spline=false)
 
void AddProperty (const G4String &key, G4MaterialPropertyVector *opv, G4bool createNewKey=false)
 
void ComputeMicroRoughnessTables ()
 
G4bool ConditionsValid (G4double E, G4double VFermi, G4double theta_i)
 
G4bool ConstPropertyExists (const char *key) const
 
G4bool ConstPropertyExists (const G4int index) const
 
G4bool ConstPropertyExists (const G4String &key) const
 
void DumpTable () const
 
 G4UCNMaterialPropertiesTable ()
 
const std::vector< std::pair< G4double, G4bool > > & GetConstProperties () const
 
G4double GetConstProperty (const char *key) const
 
G4double GetConstProperty (const G4int index) const
 
G4double GetConstProperty (const G4String &key) const
 
G4int GetConstPropertyIndex (const G4String &key) const
 
G4double GetCorrLen () const
 
const std::vector< G4String > & GetMaterialConstPropertyNames () const
 
const std::vector< G4String > & GetMaterialPropertyNames () const
 
G4doubleGetMicroRoughnessTable ()
 
G4doubleGetMicroRoughnessTransTable ()
 
G4double GetMRIntProbability (G4double, G4double)
 
G4double GetMRIntTransProbability (G4double, G4double)
 
G4double GetMRMaxProbability (G4double, G4double)
 
G4double GetMRMaxTransProbability (G4double, G4double)
 
G4double GetMRProbability (G4double, G4double, G4double, G4double, G4double)
 
G4double GetMRTransProbability (G4double, G4double, G4double, G4double, G4double)
 
const std::vector< G4MaterialPropertyVector * > & GetProperties () const
 
G4MaterialPropertyVectorGetProperty (const char *key) const
 
G4MaterialPropertyVectorGetProperty (const G4int index) const
 
G4MaterialPropertyVectorGetProperty (const G4String &key) const
 
G4int GetPropertyIndex (const G4String &key) const
 
G4double GetRMS () const
 
void InitMicroRoughnessTables ()
 
void LoadMicroRoughnessTables (G4double *, G4double *, G4double *, G4double *)
 
void RemoveConstProperty (const char *key)
 
void RemoveConstProperty (const G4String &key)
 
void RemoveProperty (const char *key)
 
void RemoveProperty (const G4String &key)
 
void SetMicroRoughnessParameters (G4double, G4double, G4int, G4int, G4double, G4double, G4double, G4double, G4int, G4int, G4double)
 
void SetMRMaxProbability (G4double, G4double, G4double)
 
void SetMRMaxTransProbability (G4double, G4double, G4double)
 
G4bool TransConditionsValid (G4double E, G4double VFermi, G4double theta_i)
 
virtual ~G4UCNMaterialPropertiesTable ()
 

Private Member Functions

G4MaterialPropertyVectorCalculateGROUPVEL ()
 

Private Attributes

G4double AngCut
 
G4double b
 
G4double E_step
 
G4double Emax
 
G4double Emin
 
std::vector< G4StringfMatConstPropNames
 
std::vector< G4StringfMatPropNames
 
std::vector< std::pair< G4double, G4bool > > fMCP
 
std::vector< G4MaterialPropertyVector * > fMP
 
G4doublemaxMicroRoughnessTable
 
G4doublemaxMicroRoughnessTransTable
 
G4int no_theta_i
 
G4int noE
 
G4doubletheMicroRoughnessTable
 
G4doubletheMicroRoughnessTransTable
 
G4double theta_i_max
 
G4double theta_i_min
 
G4double theta_i_step
 
G4double w
 

Detailed Description

Definition at line 51 of file G4UCNMaterialPropertiesTable.hh.

Constructor & Destructor Documentation

◆ G4UCNMaterialPropertiesTable()

G4UCNMaterialPropertiesTable::G4UCNMaterialPropertiesTable ( )

Definition at line 53 of file G4UCNMaterialPropertiesTable.cc.

55{
56 theMicroRoughnessTable = nullptr;
57 maxMicroRoughnessTable = nullptr;
60
62 theta_i_max = 90.*degree;
63
64 Emin = 0.e-9*eV;
65 Emax = 1000.e-9*eV;
66
67 no_theta_i = 90;
68 noE = 100;
69
71 E_step = (Emax-Emin)/(noE-1);
72
73 b = 1*nm;
74 w = 30*nm;
75
76 AngCut = 0.01*degree;
77}
static constexpr double nm
Definition: G4SIunits.hh:92
static constexpr double eV
Definition: G4SIunits.hh:201
static constexpr double degree
Definition: G4SIunits.hh:124

References AngCut, b, degree, E_step, Emax, Emin, eV, maxMicroRoughnessTable, maxMicroRoughnessTransTable, nm, no_theta_i, noE, theMicroRoughnessTable, theMicroRoughnessTransTable, theta_i_max, theta_i_min, theta_i_step, and w.

◆ ~G4UCNMaterialPropertiesTable()

G4UCNMaterialPropertiesTable::~G4UCNMaterialPropertiesTable ( )
virtual

Member Function Documentation

◆ AddConstProperty() [1/2]

void G4MaterialPropertiesTable::AddConstProperty ( const char *  key,
G4double  propertyValue,
G4bool  createNewKey = false 
)
inherited

Definition at line 446 of file G4MaterialPropertiesTable.cc.

449{
450 // Provides a way of adding a constant property to the Material Properties
451 // Table given a key
452 AddConstProperty(G4String(key), propertyValue, createNewKey);
453}
void AddConstProperty(const G4String &key, G4double propertyValue, G4bool createNewKey=false)

References G4MaterialPropertiesTable::AddConstProperty().

◆ AddConstProperty() [2/2]

void G4MaterialPropertiesTable::AddConstProperty ( const G4String key,
G4double  propertyValue,
G4bool  createNewKey = false 
)
inherited

Definition at line 417 of file G4MaterialPropertiesTable.cc.

420{
421 // Provides a way of adding a constant property to the Material Properties
422 // Table given a key
423 if(std::find(fMatConstPropNames.begin(), fMatConstPropNames.end(), key) ==
424 fMatConstPropNames.end())
425 {
426 if(createNewKey)
427 {
428 fMatConstPropNames.push_back(key);
429 fMCP.push_back(std::pair<G4double, G4bool>{ 0., true });
430 }
431 else
432 {
434 ed << "Attempting to create a new material constant property key " << key
435 << " without setting"
436 << " createNewKey parameter of AddProperty to true.";
437 G4Exception("G4MaterialPropertiesTable::AddProperty()", "mat207",
438 FatalException, ed);
439 }
440 }
441 G4int index = GetConstPropertyIndex(key);
442
443 fMCP[index] = std::pair<G4double, G4bool>{ propertyValue, true };
444}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
int G4int
Definition: G4Types.hh:85
G4int GetConstPropertyIndex(const G4String &key) const
std::vector< std::pair< G4double, G4bool > > fMCP
std::vector< G4String > fMatConstPropNames

References FatalException, G4MaterialPropertiesTable::fMatConstPropNames, G4MaterialPropertiesTable::fMCP, G4Exception(), and G4MaterialPropertiesTable::GetConstPropertyIndex().

Referenced by G4MaterialPropertiesTable::AddConstProperty(), G4GDMLReadMaterials::PropertyRead(), G4GDMLReadSolids::PropertyRead(), and SetMicroRoughnessParameters().

◆ AddEntry() [1/2]

void G4MaterialPropertiesTable::AddEntry ( const char *  key,
G4double  aPhotonEnergy,
G4double  aPropertyValue 
)
inherited

Definition at line 509 of file G4MaterialPropertiesTable.cc.

512{
513 AddEntry(G4String(key), aPhotonEnergy, aPropertyValue);
514}
void AddEntry(const G4String &key, G4double aPhotonEnergy, G4double aPropertyValue)

References G4MaterialPropertiesTable::AddEntry().

◆ AddEntry() [2/2]

void G4MaterialPropertiesTable::AddEntry ( const G4String key,
G4double  aPhotonEnergy,
G4double  aPropertyValue 
)
inherited

Definition at line 479 of file G4MaterialPropertiesTable.cc.

482{
483 // Allows to add an entry pair directly to the Material Property Vector
484 // given a key
485 if(std::find(fMatPropNames.begin(), fMatPropNames.end(), key) ==
486 fMatPropNames.end())
487 {
488 G4Exception("G4MaterialPropertiesTable::AddEntry()", "mat214",
489 FatalException, "Material Property Vector not found.");
490 }
491 G4int index = GetPropertyIndex(key);
492
493 G4MaterialPropertyVector* targetVector = fMP[index];
494 if(targetVector != nullptr)
495 {
496 targetVector->InsertValues(aPhotonEnergy, aPropertyValue);
497 }
498 else
499 {
500 G4Exception("G4MaterialPropertiesTable::AddEntry()", "mat208",
501 FatalException, "Material Property Vector not found.");
502 }
503 if(key == "RINDEX")
504 {
506 }
507}
G4MaterialPropertyVector * CalculateGROUPVEL()
G4int GetPropertyIndex(const G4String &key) const
std::vector< G4MaterialPropertyVector * > fMP
std::vector< G4String > fMatPropNames
void InsertValues(const G4double energy, const G4double value)

References G4MaterialPropertiesTable::CalculateGROUPVEL(), FatalException, G4MaterialPropertiesTable::fMatPropNames, G4MaterialPropertiesTable::fMP, G4Exception(), G4MaterialPropertiesTable::GetPropertyIndex(), and G4PhysicsFreeVector::InsertValues().

Referenced by G4MaterialPropertiesTable::AddEntry().

◆ AddProperty() [1/5]

G4MaterialPropertyVector * G4MaterialPropertiesTable::AddProperty ( const char *  key,
G4double photonEnergies,
G4double propertyValues,
G4int  numEntries,
G4bool  createNewKey = false,
G4bool  spline = false 
)
inherited

Definition at line 351 of file G4MaterialPropertiesTable.cc.

354{
355
356 // Provides a way of adding a property to the Material Properties
357 // Table given a pair of arrays and a key
358 G4String k(key);
359
360 std::vector<G4double> energies(photonEnergies, photonEnergies + numEntries);
361 std::vector<G4double> values(propertyValues, propertyValues + numEntries);
362 return AddProperty(k, energies, values, createNewKey, spline);
363}
G4MaterialPropertyVector * AddProperty(const G4String &key, const std::vector< G4double > &photonEnergies, const std::vector< G4double > &propertyValues, G4bool createNewKey=false, G4bool spline=false)

References G4MaterialPropertiesTable::AddProperty().

◆ AddProperty() [2/5]

void G4MaterialPropertiesTable::AddProperty ( const char *  key,
G4MaterialPropertyVector opv,
G4bool  createNewKey = false 
)
inherited

Definition at line 401 of file G4MaterialPropertiesTable.cc.

404{
405 AddProperty(G4String(key), mpv, createNewKey);
406}

References G4MaterialPropertiesTable::AddProperty().

◆ AddProperty() [3/5]

void G4MaterialPropertiesTable::AddProperty ( const G4String key,
const G4String mat 
)
inherited

Definition at line 408 of file G4MaterialPropertiesTable.cc.

410{
411 // load a material property vector defined in Geant4 source
414 AddProperty(key, v);
415}
G4MaterialPropertyVector * GetProperty(const G4String &key, const G4String &mat)

References G4MaterialPropertiesTable::AddProperty(), and G4OpticalMaterialProperties::GetProperty().

◆ AddProperty() [4/5]

G4MaterialPropertyVector * G4MaterialPropertiesTable::AddProperty ( const G4String key,
const std::vector< G4double > &  photonEnergies,
const std::vector< G4double > &  propertyValues,
G4bool  createNewKey = false,
G4bool  spline = false 
)
inherited

Definition at line 298 of file G4MaterialPropertiesTable.cc.

302{
303 if(photonEnergies.size() != propertyValues.size())
304 {
306 ed << "AddProperty error!";
307 G4Exception("G4MaterialPropertiesTable::AddProperty()", "mat204",
308 FatalException, ed);
309 }
310
311 // if the key doesn't exist, add it if requested
312 if(std::find(fMatPropNames.begin(), fMatPropNames.end(), key) ==
313 fMatPropNames.end())
314 {
315 if(createNewKey)
316 {
317 fMatPropNames.push_back(key);
318 fMP.push_back(nullptr);
319 }
320 else
321 {
323 ed << "Attempting to create a new material property key " << key
324 << " without setting\n"
325 << "createNewKey parameter of AddProperty to true.";
326 G4Exception("G4MaterialPropertiesTable::AddProperty()", "mat205",
327 FatalException, ed);
328 }
329 }
330
332 new G4MaterialPropertyVector(photonEnergies, propertyValues, spline);
333 mpv->SetVerboseLevel(1);
334 if(spline)
335 {
337 }
338 G4int index = GetPropertyIndex(key);
339 fMP[index] = mpv;
340
341 // if key is RINDEX, we calculate GROUPVEL -
342 // contribution from Tao Lin (IHEP, the JUNO experiment)
343 if(key == "RINDEX")
344 {
346 }
347
348 return mpv;
349}
G4PhysicsFreeVector G4MaterialPropertyVector
void SetVerboseLevel(G4int value)
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)

References G4MaterialPropertiesTable::CalculateGROUPVEL(), FatalException, G4PhysicsVector::FillSecondDerivatives(), G4MaterialPropertiesTable::fMatPropNames, G4MaterialPropertiesTable::fMP, G4Exception(), G4MaterialPropertiesTable::GetPropertyIndex(), and G4PhysicsVector::SetVerboseLevel().

Referenced by G4MaterialPropertiesTable::AddProperty(), G4MaterialPropertiesTable::CalculateGROUPVEL(), G4GDMLReadMaterials::PropertyRead(), and G4GDMLReadSolids::PropertyRead().

◆ AddProperty() [5/5]

void G4MaterialPropertiesTable::AddProperty ( const G4String key,
G4MaterialPropertyVector opv,
G4bool  createNewKey = false 
)
inherited

Definition at line 365 of file G4MaterialPropertiesTable.cc.

368{
369 // Provides a way of adding a property to the Material Properties
370 // Table given an G4MaterialPropertyVector Reference and a key
371 // if the key doesn't exist, add it
372 if(std::find(fMatPropNames.begin(), fMatPropNames.end(), key) ==
373 fMatPropNames.end())
374 {
375 if(createNewKey)
376 {
377 fMatPropNames.push_back(key);
378 fMP.push_back(nullptr);
379 }
380 else
381 {
383 ed << "Attempting to create a new material property key " << key
384 << " without setting\n"
385 << "createNewKey parameter of AddProperty to true.";
386 G4Exception("G4MaterialPropertiesTable::AddProperty()", "mat206",
387 FatalException, ed);
388 }
389 }
390 G4int index = GetPropertyIndex(key);
391 fMP[index] = mpv;
392
393 // if key is RINDEX, we calculate GROUPVEL -
394 // contribution from Tao Lin (IHEP, the JUNO experiment)
395 if(key == "RINDEX")
396 {
398 }
399}

References G4MaterialPropertiesTable::CalculateGROUPVEL(), FatalException, G4MaterialPropertiesTable::fMatPropNames, G4MaterialPropertiesTable::fMP, G4Exception(), and G4MaterialPropertiesTable::GetPropertyIndex().

◆ CalculateGROUPVEL()

G4MaterialPropertyVector * G4MaterialPropertiesTable::CalculateGROUPVEL ( )
privateinherited

Definition at line 542 of file G4MaterialPropertiesTable.cc.

543{
544#ifdef G4MULTITHREADED
545 G4AutoLock mptm(&materialPropertyTableMutex);
546#endif
547
548 // check if "GROUPVEL" already exists. If so, remove it.
549 if(fMP[kGROUPVEL] != nullptr)
550 {
551 this->RemoveProperty("GROUPVEL");
552 }
553
554 // fetch RINDEX data, give up if unavailable
556 if(rindex == nullptr)
557 {
558 return nullptr;
559 }
560
561 // RINDEX exists but has no entries, give up
562 if(rindex->GetVectorLength() == 0)
563 {
564 return nullptr;
565 }
566
567 // add GROUPVEL vector
569 groupvel->SetVerboseLevel(1);
570
571 // fill GROUPVEL vector using RINDEX values
572 // rindex built-in "iterator" was advanced to first entry above
573 G4double E0 = rindex->Energy(0);
574 G4double n0 = (*rindex)[0];
575
576 if(E0 <= 0.)
577 {
578 G4Exception("G4MaterialPropertiesTable::CalculateGROUPVEL()", "mat211",
579 FatalException, "Optical Photon Energy <= 0");
580 }
581
582 if(rindex->GetVectorLength() >= 2)
583 {
584 // good, we have at least two entries in RINDEX
585 // get next energy/value pair
586
587 G4double E1 = rindex->Energy(1);
588 G4double n1 = (*rindex)[1];
589
590 if(E1 <= 0.)
591 {
592 G4Exception("G4MaterialPropertiesTable::CalculateGROUPVEL()", "mat212",
593 FatalException, "Optical Photon Energy <= 0");
594 }
595
596 G4double vg;
597
598 // add entry at first photon energy
599 vg = c_light / (n0 + (n1 - n0) / G4Log(E1 / E0));
600
601 // allow only for 'normal dispersion' -> dn/d(logE) > 0
602 if((vg < 0) || (vg > c_light / n0))
603 {
604 vg = c_light / n0;
605 }
606
607 groupvel->InsertValues(E0, vg);
608
609 // add entries at midpoints between remaining photon energies
610 for(size_t i = 2; i < rindex->GetVectorLength(); ++i)
611 {
612 vg = c_light / (0.5 * (n0 + n1) + (n1 - n0) / G4Log(E1 / E0));
613
614 // allow only for 'normal dispersion' -> dn/d(logE) > 0
615 if((vg < 0) || (vg > c_light / (0.5 * (n0 + n1))))
616 {
617 vg = c_light / (0.5 * (n0 + n1));
618 }
619 groupvel->InsertValues(0.5 * (E0 + E1), vg);
620
621 // get next energy/value pair, or exit loop
622 E0 = E1;
623 n0 = n1;
624 E1 = rindex->Energy(i);
625 n1 = (*rindex)[i];
626
627 if(E1 <= 0.)
628 {
629 G4Exception("G4MaterialPropertiesTable::CalculateGROUPVEL()", "mat213",
630 FatalException, "Optical Photon Energy <= 0");
631 }
632 }
633
634 // add entry at last photon energy
635 vg = c_light / (n1 + (n1 - n0) / G4Log(E1 / E0));
636
637 // allow only for 'normal dispersion' -> dn/d(logE) > 0
638 if((vg < 0) || (vg > c_light / n1))
639 {
640 vg = c_light / n1;
641 }
642 groupvel->InsertValues(E1, vg);
643 }
644 else // only one entry in RINDEX -- weird!
645 {
646 groupvel->InsertValues(E0, c_light / n0);
647 }
648
649 this->AddProperty("GROUPVEL", groupvel);
650
651 return groupvel;
652}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83
void RemoveProperty(const G4String &key)
G4MaterialPropertyVector * GetProperty(const char *key) const
G4double Energy(const std::size_t index) const
std::size_t GetVectorLength() const
float c_light
Definition: hepunit.py:256

References G4MaterialPropertiesTable::AddProperty(), source.hepunit::c_light, G4PhysicsVector::Energy(), FatalException, G4MaterialPropertiesTable::fMP, G4Exception(), G4Log(), G4MaterialPropertiesTable::GetProperty(), G4PhysicsVector::GetVectorLength(), G4PhysicsFreeVector::InsertValues(), kGROUPVEL, kRINDEX, G4MaterialPropertiesTable::RemoveProperty(), and G4PhysicsVector::SetVerboseLevel().

Referenced by G4MaterialPropertiesTable::AddEntry(), and G4MaterialPropertiesTable::AddProperty().

◆ ComputeMicroRoughnessTables()

void G4UCNMaterialPropertiesTable::ComputeMicroRoughnessTables ( )

Definition at line 141 of file G4UCNMaterialPropertiesTable.cc.

142{
143// Reads the parameters for the mr-probability computation from the
144// corresponding material properties and stores it in the appropriate
145// variables
146
147 b = GetConstProperty("MR_RRMS");
148 G4double b2 = b*b;
149 w = GetConstProperty("MR_CORRLEN");
150 G4double w2 = w*w;
151
152 no_theta_i = G4int(GetConstProperty("MR_NBTHETA")+0.1);
153 //G4cout << "no_theta: " << no_theta_i << G4endl;
154 noE = G4int(GetConstProperty("MR_NBE")+0.1);
155 //G4cout << "noE: " << noE << G4endl;
156
157 theta_i_min = GetConstProperty("MR_THETAMIN");
158 theta_i_max = GetConstProperty("MR_THETAMAX");
159 Emin = GetConstProperty("MR_EMIN");
160 Emax = GetConstProperty("MR_EMAX");
161 G4int AngNoTheta = G4int(GetConstProperty("MR_ANGNOTHETA")+0.1);
162 G4int AngNoPhi = G4int(GetConstProperty("MR_ANGNOPHI")+0.1);
163 AngCut = GetConstProperty("MR_ANGCUT");
164
165 // The Fermi potential was saved in neV by P.F.
166
167 G4double fermipot = GetConstProperty("FERMIPOT")*(1.e-9*eV);
168
169 //G4cout << "Fermipot: " << fermipot/(1.e-9*eV) << "neV" << G4endl;
170
171 G4double theta_i, E;
172
173 // Calculates the increment in theta_i in the lookup-table
175
176 //G4cout << "theta_i_step: " << theta_i_step << G4endl;
177
178 // Calculates the increment in energy in the lookup-table
179 E_step = (Emax-Emin)/(noE-1);
180
181 // Runs the lookup-table memory allocation
183
184 G4int counter = 0;
185
186 //G4cout << "Calculating Microroughnesstables...";
187
188 // Writes the mr-lookup-tables to files for immediate control
189
190 std::ofstream dateir("MRrefl.dat",std::ios::out);
191 std::ofstream dateit("MRtrans.dat",std::ios::out);
192
193 //G4cout << theMicroRoughnessTable << G4endl;
194
195 for (theta_i=theta_i_min; theta_i<=theta_i_max+1e-6; theta_i+=theta_i_step) {
196 // Calculation for each cell in the lookup-table
197 for (E=Emin; E<=Emax; E+=E_step) {
198 *(theMicroRoughnessTable+counter) =
200 IntIplus(E, fermipot, theta_i, AngNoTheta, AngNoPhi,
201 b2, w2, maxMicroRoughnessTable+counter, AngCut);
202
203 *(theMicroRoughnessTransTable+counter) =
205 IntIminus(E, fermipot, theta_i, AngNoTheta, AngNoPhi,
206 b2, w2, maxMicroRoughnessTransTable+counter,
207 AngCut);
208
209 dateir << *(theMicroRoughnessTable+counter) << G4endl;
210 dateit << *(theMicroRoughnessTransTable+counter) << G4endl;
211
212 counter++;
213
214 //G4cout << counter << G4endl;
215 }
216 }
217
218 dateir.close();
219 dateit.close();
220
221 // Writes the mr-lookup-tables to files for immediate control
222
223 std::ofstream dateic("MRcheck.dat",std::ios::out);
224 std::ofstream dateimr("MRmaxrefl.dat",std::ios::out);
225 std::ofstream dateimt("MRmaxtrans.dat",std::ios::out);
226
227 for (theta_i=theta_i_min; theta_i<=theta_i_max+1e-6; theta_i+=theta_i_step) {
228 for (E=Emin; E<=Emax; E+=E_step) {
229
230 // tests the GetXXProbability functions by writing the entries
231 // of the lookup tables to files
232
233 dateic << GetMRIntProbability(theta_i,E) << G4endl;
234 dateimr << GetMRMaxProbability(theta_i,E) << G4endl;
235 dateimt << GetMRMaxTransProbability(theta_i,E) << G4endl;
236 }
237 }
238
239 dateic.close();
240 dateimr.close();
241 dateimt.close();
242}
#define G4endl
Definition: G4ios.hh:57
G4double GetConstProperty(const G4String &key) const
G4double GetMRMaxTransProbability(G4double, G4double)
G4double GetMRMaxProbability(G4double, G4double)
G4double GetMRIntProbability(G4double, G4double)
static G4UCNMicroRoughnessHelper * GetInstance()

References AngCut, b, E_step, Emax, Emin, eV, G4endl, G4MaterialPropertiesTable::GetConstProperty(), G4UCNMicroRoughnessHelper::GetInstance(), GetMRIntProbability(), GetMRMaxProbability(), GetMRMaxTransProbability(), InitMicroRoughnessTables(), maxMicroRoughnessTable, maxMicroRoughnessTransTable, no_theta_i, noE, theMicroRoughnessTable, theMicroRoughnessTransTable, theta_i_max, theta_i_min, theta_i_step, and w.

Referenced by SetMicroRoughnessParameters().

◆ ConditionsValid()

G4bool G4UCNMaterialPropertiesTable::ConditionsValid ( G4double  E,
G4double  VFermi,
G4double  theta_i 
)

Definition at line 413 of file G4UCNMaterialPropertiesTable.cc.

416{
417 G4double k = std::sqrt(2*neutron_mass_c2*E / hbarc_squared);
418 G4double k_l = std::sqrt(2*neutron_mass_c2*VFermi / hbarc_squared);
419
420 //G4cout << " Energy: " << E/(1.e-9*eV) << "neV"
421 // << " VFermi: " << VFermi/(1.e-9*eV) << "neV"
422 // << " theta: " << theta_i/degree << "degree" << G4endl;
423
424 //G4cout << " Conditions - 2*b*k*cos(theta): " << 2*b*k*cos(theta_i)
425 // << ", 2*b*kl: " << 2*b*k_l << G4endl;
426
427 // see eq. 17 of the Steyerl paper
428
429 if (2*b*k*std::cos(theta_i) < 1 && 2*b*k_l < 1) return true;
430 else return false;
431}
float hbarc_squared
Definition: hepunit.py:265
float neutron_mass_c2
Definition: hepunit.py:275

References b, source.hepunit::hbarc_squared, and source.hepunit::neutron_mass_c2.

◆ ConstPropertyExists() [1/3]

G4bool G4MaterialPropertiesTable::ConstPropertyExists ( const char *  key) const
inherited

Definition at line 253 of file G4MaterialPropertiesTable.cc.

254{
255 size_t index = std::distance(
256 fMatConstPropNames.begin(),
257 std::find(fMatConstPropNames.begin(), fMatConstPropNames.end(), key));
258 if(index < fMatConstPropNames.size()) // index is type size_t so >= 0
259 return ConstPropertyExists(index);
260 return false;
261}
G4bool ConstPropertyExists(const G4String &key) const

References G4MaterialPropertiesTable::ConstPropertyExists(), and G4MaterialPropertiesTable::fMatConstPropNames.

◆ ConstPropertyExists() [2/3]

G4bool G4MaterialPropertiesTable::ConstPropertyExists ( const G4int  index) const
inherited

Definition at line 231 of file G4MaterialPropertiesTable.cc.

232{
233 // Returns true if a const property corresponding to 'index' exists
234
235 if(index >= 0 && index < (G4int) fMCP.size() && fMCP[index].second == true)
236 {
237 return true;
238 }
239 return false;
240}

References G4MaterialPropertiesTable::fMCP.

◆ ConstPropertyExists() [3/3]

G4bool G4MaterialPropertiesTable::ConstPropertyExists ( const G4String key) const
inherited

◆ DumpTable()

void G4MaterialPropertiesTable::DumpTable ( ) const
inherited

Definition at line 516 of file G4MaterialPropertiesTable.cc.

517{
518 // material properties
519 G4int j = 0;
520 for(const auto& prop : fMP)
521 {
522 if(prop != nullptr)
523 {
524 G4cout << j << ": " << fMatPropNames[j] << G4endl;
525 prop->DumpValues();
526 }
527 ++j;
528 }
529 // material constant properties
530 j = 0;
531 for(const auto& cprop : fMCP)
532 {
533 if(cprop.second == true)
534 {
535 G4cout << j << ": " << fMatConstPropNames[j] << " " << cprop.first
536 << G4endl;
537 }
538 ++j;
539 }
540}
G4GLOB_DLL std::ostream G4cout

References G4MaterialPropertiesTable::fMatConstPropNames, G4MaterialPropertiesTable::fMatPropNames, G4MaterialPropertiesTable::fMCP, G4MaterialPropertiesTable::fMP, G4cout, and G4endl.

◆ GetConstProperties()

const std::vector< std::pair< G4double, G4bool > > & G4MaterialPropertiesTable::GetConstProperties ( ) const
inlineinherited

◆ GetConstProperty() [1/3]

G4double G4MaterialPropertiesTable::GetConstProperty ( const char *  key) const
inherited

◆ GetConstProperty() [2/3]

G4double G4MaterialPropertiesTable::GetConstProperty ( const G4int  index) const
inherited

Definition at line 204 of file G4MaterialPropertiesTable.cc.

205{
206 // Returns the constant material property corresponding to an index
207 // fatal exception if property not found
208
209 if(index < (G4int) fMCP.size() && fMCP[index].second == true)
210 return fMCP[index].first;
212 ed << "Constant Material Property Index " << index << " not found.";
213 G4Exception("G4MaterialPropertiesTable::GetConstProperty()", "mat202",
214 FatalException, ed);
215 return 0.;
216}

References FatalException, G4MaterialPropertiesTable::fMCP, and G4Exception().

◆ GetConstProperty() [3/3]

G4double G4MaterialPropertiesTable::GetConstProperty ( const G4String key) const
inherited

◆ GetConstPropertyIndex()

G4int G4MaterialPropertiesTable::GetConstPropertyIndex ( const G4String key) const
inherited

Definition at line 171 of file G4MaterialPropertiesTable.cc.

173{
174 // Returns the constant material property index corresponding to a key
175
176 size_t index = std::distance(
177 fMatConstPropNames.begin(),
178 std::find(fMatConstPropNames.begin(), fMatConstPropNames.end(), key));
179 if(index < fMatConstPropNames.size())
180 return index;
181
183 ed << "Constant Material Property Index for key " << key << " not found.";
184 G4Exception("G4MaterialPropertiesTable::GetConstPropertyIndex()", "mat200",
185 FatalException, ed);
186 return 0;
187}

References FatalException, G4MaterialPropertiesTable::fMatConstPropNames, and G4Exception().

Referenced by G4MaterialPropertiesTable::AddConstProperty(), G4MaterialPropertiesTable::GetConstProperty(), and G4MaterialPropertiesTable::RemoveConstProperty().

◆ GetCorrLen()

G4double G4UCNMaterialPropertiesTable::GetCorrLen ( ) const
inline

Definition at line 181 of file G4UCNMaterialPropertiesTable.hh.

181{return w;}

References w.

Referenced by G4UCNBoundaryProcess::Loss().

◆ GetMaterialConstPropertyNames()

const std::vector< G4String > & G4MaterialPropertiesTable::GetMaterialConstPropertyNames ( ) const
inlineinherited

◆ GetMaterialPropertyNames()

const std::vector< G4String > & G4MaterialPropertiesTable::GetMaterialPropertyNames ( ) const
inlineinherited

◆ GetMicroRoughnessTable()

G4double * G4UCNMaterialPropertiesTable::GetMicroRoughnessTable ( )

Definition at line 87 of file G4UCNMaterialPropertiesTable.cc.

88{
90}

References theMicroRoughnessTable.

Referenced by G4UCNBoundaryProcess::PostStepDoIt().

◆ GetMicroRoughnessTransTable()

G4double * G4UCNMaterialPropertiesTable::GetMicroRoughnessTransTable ( )

Definition at line 92 of file G4UCNMaterialPropertiesTable.cc.

93{
95}

References theMicroRoughnessTransTable.

◆ GetMRIntProbability()

G4double G4UCNMaterialPropertiesTable::GetMRIntProbability ( G4double  theta_i,
G4double  Energy 
)

Definition at line 244 of file G4UCNMaterialPropertiesTable.cc.

246{
248 G4cout << "Do not have theMicroRoughnessTable" << G4endl;
249 return 0.;
250 }
251
252 // if theta_i or energy outside the range for which the lookup table is
253 // calculated, the probability is set to zero
254
255 //G4cout << "theta_i: " << theta_i/degree << "degree"
256 // << " theta_i_min: " << theta_i_min/degree << "degree"
257 // << " theta_i_max: " << theta_i_max/degree << "degree"
258 // << " Energy: " << Energy/(1.e-9*eV) << "neV"
259 // << " Emin: " << Emin/(1.e-9*eV) << "neV"
260 // << " Emax: " << Emax/(1.e-9*eV) << "neV" << G4endl;
261
262 if (theta_i<theta_i_min || theta_i>theta_i_max || Energy<Emin || Energy>Emax)
263 return 0.;
264
265 // Determines the nearest cell in the lookup table which contains
266 // the probability
267
268 G4int theta_i_pos = G4int((theta_i-theta_i_min)/theta_i_step+0.5);
269 G4int E_pos = G4int((Energy-Emin)/E_step+0.5);
270
271 // lookup table is onedimensional (1 row), energy is in rows,
272 // theta_i in columns
273
274 //G4cout << "E_pos: " << E_pos << " theta_i_pos: " << theta_i_pos << G4endl;
275 //G4cout << "Probability: " << *(theMicroRoughnessTable+E_pos+theta_i_pos*(noE-1)) << G4endl;
276
277 return *(theMicroRoughnessTable+E_pos+theta_i_pos*(noE - 1));
278}

References E_step, Emax, Emin, G4cout, G4endl, noE, theMicroRoughnessTable, theta_i_max, theta_i_min, and theta_i_step.

Referenced by ComputeMicroRoughnessTables().

◆ GetMRIntTransProbability()

G4double G4UCNMaterialPropertiesTable::GetMRIntTransProbability ( G4double  theta_i,
G4double  Energy 
)

Definition at line 280 of file G4UCNMaterialPropertiesTable.cc.

282{
283 if (!theMicroRoughnessTransTable) return 0.;
284
285 // if theta_i or energy outside the range for which the lookup table
286 // is calculated, the probability is set to zero
287
288 if (theta_i<theta_i_min || theta_i>theta_i_max || Energy<Emin || Energy>Emax)
289 return 0.;
290
291 // Determines the nearest cell in the lookup table which contains
292 // the probability
293
294 G4int theta_i_pos = G4int((theta_i-theta_i_min)/theta_i_step+0.5);
295 G4int E_pos = G4int((Energy-Emin)/E_step+0.5);
296
297 // lookup table is onedimensional (1 row), energy is in rows,
298 // theta_i in columns
299
300 return *(theMicroRoughnessTransTable+E_pos+theta_i_pos*(noE - 1));
301}

References E_step, Emax, Emin, noE, theMicroRoughnessTransTable, theta_i_max, theta_i_min, and theta_i_step.

◆ GetMRMaxProbability()

G4double G4UCNMaterialPropertiesTable::GetMRMaxProbability ( G4double  theta_i,
G4double  Energy 
)

Definition at line 303 of file G4UCNMaterialPropertiesTable.cc.

305{
306 if (!maxMicroRoughnessTable) return 0.;
307
308 // if theta_i or energy outside the range for which the lookup table
309 // is calculated, the probability is set to zero
310
311 if (theta_i<theta_i_min || theta_i>theta_i_max || Energy<Emin || Energy>Emax)
312 return 0.;
313
314 // Determines the nearest cell in the lookup table which contains
315 // the probability
316
317 G4int theta_i_pos = G4int((theta_i-theta_i_min)/theta_i_step+0.5);
318 G4int E_pos = G4int((Energy-Emin)/E_step+0.5);
319
320 // lookup table is onedimensional (1 row), energy is in rows,
321 // theta_i in columns
322
323 return *(maxMicroRoughnessTable+E_pos+theta_i_pos*noE);
324}

References E_step, Emax, Emin, maxMicroRoughnessTable, noE, theta_i_max, theta_i_min, and theta_i_step.

Referenced by ComputeMicroRoughnessTables().

◆ GetMRMaxTransProbability()

G4double G4UCNMaterialPropertiesTable::GetMRMaxTransProbability ( G4double  theta_i,
G4double  Energy 
)

Definition at line 349 of file G4UCNMaterialPropertiesTable.cc.

351{
352 if (!maxMicroRoughnessTransTable) return 0.;
353
354 // if theta_i or energy outside the range for which the lookup table
355 // is calculated, the probability is set to zero
356
357 if (theta_i<theta_i_min || theta_i>theta_i_max || Energy<Emin || Energy>Emax)
358 return 0.;
359
360 // Determines the nearest cell in the lookup table which contains
361 // the probability
362
363 G4int theta_i_pos = G4int((theta_i-theta_i_min)/theta_i_step+0.5);
364 G4int E_pos = G4int((Energy-Emin)/E_step+0.5);
365
366 // lookup table is onedimensional (1 row), energy is in rows,
367 // theta_i in columns
368
369 return *(maxMicroRoughnessTransTable+E_pos+theta_i_pos*noE);
370}

References E_step, Emax, Emin, maxMicroRoughnessTransTable, noE, theta_i_max, theta_i_min, and theta_i_step.

Referenced by ComputeMicroRoughnessTables().

◆ GetMRProbability()

G4double G4UCNMaterialPropertiesTable::GetMRProbability ( G4double  theta_i,
G4double  Energy,
G4double  fermipot,
G4double  theta_o,
G4double  phi_o 
)

Definition at line 395 of file G4UCNMaterialPropertiesTable.cc.

399{
401 ProbIplus(Energy, fermipot, theta_i, theta_o, phi_o, b, w, AngCut);
402}

References AngCut, b, G4UCNMicroRoughnessHelper::GetInstance(), and w.

◆ GetMRTransProbability()

G4double G4UCNMaterialPropertiesTable::GetMRTransProbability ( G4double  theta_i,
G4double  Energy,
G4double  fermipot,
G4double  theta_o,
G4double  phi_o 
)

Definition at line 404 of file G4UCNMaterialPropertiesTable.cc.

408{
410 ProbIminus(Energy, fermipot,theta_i, theta_o, phi_o, b, w, AngCut);
411}

References AngCut, b, G4UCNMicroRoughnessHelper::GetInstance(), and w.

◆ GetProperties()

const std::vector< G4MaterialPropertyVector * > & G4MaterialPropertiesTable::GetProperties ( ) const
inlineinherited

Definition at line 155 of file G4MaterialPropertiesTable.hh.

156 {
157 return fMP;
158 }

References G4MaterialPropertiesTable::fMP.

Referenced by G4GDMLWriteMaterials::PropertyWrite(), and G4GDMLWriteSolids::PropertyWrite().

◆ GetProperty() [1/3]

G4MaterialPropertyVector * G4MaterialPropertiesTable::GetProperty ( const char *  key) const
inherited

◆ GetProperty() [2/3]

G4MaterialPropertyVector * G4MaterialPropertiesTable::GetProperty ( const G4int  index) const
inherited

Definition at line 288 of file G4MaterialPropertiesTable.cc.

290{
291 // Returns a Material Property Vector corresponding to an index
292 // returns nullptr if the property has not been defined by user
293 if(index >= 0 && index < (G4int) fMP.size())
294 return fMP[index];
295 return nullptr;
296}

References G4MaterialPropertiesTable::fMP.

◆ GetProperty() [3/3]

G4MaterialPropertyVector * G4MaterialPropertiesTable::GetProperty ( const G4String key) const
inherited

Definition at line 263 of file G4MaterialPropertiesTable.cc.

265{
266 // Returns a Material Property Vector corresponding to a key
267 if(std::find(fMatPropNames.begin(), fMatPropNames.end(), key) !=
268 fMatPropNames.end())
269 {
270 const G4int index = GetPropertyIndex(G4String(key));
271 return GetProperty(index);
272 }
273 return nullptr;
274}

References G4MaterialPropertiesTable::fMatPropNames, G4MaterialPropertiesTable::GetProperty(), and G4MaterialPropertiesTable::GetPropertyIndex().

◆ GetPropertyIndex()

G4int G4MaterialPropertiesTable::GetPropertyIndex ( const G4String key) const
inherited

Definition at line 189 of file G4MaterialPropertiesTable.cc.

190{
191 // Returns the material property index corresponding to a key
192 size_t index =
193 std::distance(fMatPropNames.begin(),
194 std::find(fMatPropNames.begin(), fMatPropNames.end(), key));
195 if(index < fMatPropNames.size())
196 return index;
198 ed << "Material Property Index for key " << key << " not found.";
199 G4Exception("G4MaterialPropertiesTable::GetPropertyIndex()", "mat201",
200 FatalException, ed);
201 return 0;
202}

References FatalException, G4MaterialPropertiesTable::fMatPropNames, and G4Exception().

Referenced by G4MaterialPropertiesTable::AddEntry(), G4MaterialPropertiesTable::AddProperty(), G4MaterialPropertiesTable::GetProperty(), and G4MaterialPropertiesTable::RemoveProperty().

◆ GetRMS()

G4double G4UCNMaterialPropertiesTable::GetRMS ( ) const
inline

Definition at line 180 of file G4UCNMaterialPropertiesTable.hh.

180{return b;}

References b.

Referenced by G4UCNBoundaryProcess::Loss().

◆ InitMicroRoughnessTables()

void G4UCNMaterialPropertiesTable::InitMicroRoughnessTables ( )

Definition at line 109 of file G4UCNMaterialPropertiesTable.cc.

110{
111 G4int NEdim = 0;
112 G4int Nthetadim = 0;
113
114 // Checks if the number of angles is available and stores it
115
116 if (ConstPropertyExists("MR_NBTHETA"))
117 Nthetadim = G4int(GetConstProperty("MR_NBTHETA")+0.1);
118
119 // Checks if the number of energies is available and stores it
120
121 if (ConstPropertyExists("MR_NBE"))
122 NEdim = G4int(GetConstProperty("MR_NBE")+0.1);
123
124 //G4cout << "thetadim: " << Nthetadim << " , Edim: " << NEdim << G4endl;
125
126 // If both dimensions of the lookup-table are non-trivial:
127 // delete old tables if existing and allocate memory for new tables
128
129 if (Nthetadim*NEdim > 0) {
131 theMicroRoughnessTable = new G4double[Nthetadim*NEdim];
133 maxMicroRoughnessTable = new G4double[Nthetadim*NEdim];
135 theMicroRoughnessTransTable = new G4double[Nthetadim*NEdim];
137 maxMicroRoughnessTransTable = new G4double[Nthetadim*NEdim];
138 }
139}

References G4MaterialPropertiesTable::ConstPropertyExists(), G4MaterialPropertiesTable::GetConstProperty(), maxMicroRoughnessTable, maxMicroRoughnessTransTable, theMicroRoughnessTable, and theMicroRoughnessTransTable.

Referenced by ComputeMicroRoughnessTables().

◆ LoadMicroRoughnessTables()

void G4UCNMaterialPropertiesTable::LoadMicroRoughnessTables ( G4double pMicroRoughnessTable,
G4double pmaxMicroRoughnessTable,
G4double pMicroRoughnessTransTable,
G4double pmaxMicroRoughnessTransTable 
)

Definition at line 97 of file G4UCNMaterialPropertiesTable.cc.

102{
103 theMicroRoughnessTable = pMicroRoughnessTable;
104 maxMicroRoughnessTable = pmaxMicroRoughnessTable;
105 theMicroRoughnessTransTable = pMicroRoughnessTransTable;
106 maxMicroRoughnessTransTable = pmaxMicroRoughnessTransTable;
107}

References maxMicroRoughnessTable, maxMicroRoughnessTransTable, theMicroRoughnessTable, and theMicroRoughnessTransTable.

◆ RemoveConstProperty() [1/2]

void G4MaterialPropertiesTable::RemoveConstProperty ( const char *  key)
inherited

Definition at line 462 of file G4MaterialPropertiesTable.cc.

463{
465}
void RemoveConstProperty(const G4String &key)

References G4MaterialPropertiesTable::RemoveConstProperty().

◆ RemoveConstProperty() [2/2]

void G4MaterialPropertiesTable::RemoveConstProperty ( const G4String key)
inherited

Definition at line 455 of file G4MaterialPropertiesTable.cc.

456{
457 G4int index = GetConstPropertyIndex(key);
458 if(index < (G4int) fMCP.size())
459 fMCP[index] = std::pair<G4double, G4bool>{ 0., false };
460}

References G4MaterialPropertiesTable::fMCP, and G4MaterialPropertiesTable::GetConstPropertyIndex().

Referenced by G4MaterialPropertiesTable::RemoveConstProperty(), and SetMicroRoughnessParameters().

◆ RemoveProperty() [1/2]

void G4MaterialPropertiesTable::RemoveProperty ( const char *  key)
inherited

Definition at line 474 of file G4MaterialPropertiesTable.cc.

475{
477}

References G4MaterialPropertiesTable::RemoveProperty().

◆ RemoveProperty() [2/2]

void G4MaterialPropertiesTable::RemoveProperty ( const G4String key)
inherited

◆ SetMicroRoughnessParameters()

void G4UCNMaterialPropertiesTable::SetMicroRoughnessParameters ( G4double  ww,
G4double  bb,
G4int  no_theta,
G4int  no_E,
G4double  theta_min,
G4double  theta_max,
G4double  E_min,
G4double  E_max,
G4int  AngNoTheta,
G4int  AngNoPhi,
G4double  AngularCut 
)

Definition at line 453 of file G4UCNMaterialPropertiesTable.cc.

460{
461 //G4cout << "Setting Microroughness Parameters...";
462
463 // Removes an existing RMS roughness
464 if (ConstPropertyExists("MR_RRMS")) RemoveConstProperty("MR_RRMS");
465
466 // Adds a new RMS roughness
467 AddConstProperty("MR_RRMS", bb);
468
469 //G4cout << "b: " << bb << G4endl;
470
471 // Removes an existing correlation length
472 if (ConstPropertyExists("MR_CORRLEN")) RemoveConstProperty("MR_CORRLEN");
473
474 // Adds a new correlation length
475 AddConstProperty("MR_CORRLEN", ww);
476
477 //G4cout << "w: " << ww << G4endl;
478
479 // Removes an existing number of thetas
480 if (ConstPropertyExists("MR_NBTHETA")) RemoveConstProperty("MR_NBTHETA");
481
482 // Adds a new number of thetas
483 AddConstProperty("MR_NBTHETA", (G4double)no_theta);
484
485 //G4cout << "no_theta: " << no_theta << G4endl;
486
487 // Removes an existing number of Energies
488 if (ConstPropertyExists("MR_NBE")) RemoveConstProperty("MR_NBE");
489
490 // Adds a new number of Energies
491 AddConstProperty("MR_NBE", (G4double)no_E);
492
493 //G4cout << "no_E: " << no_E << G4endl;
494
495 // Removes an existing minimum theta
496 if (ConstPropertyExists("MR_THETAMIN")) RemoveConstProperty("MR_THETAMIN");
497
498 // Adds a new minimum theta
499 AddConstProperty("MR_THETAMIN", theta_min);
500
501 //G4cout << "theta_min: " << theta_min << G4endl;
502
503 // Removes an existing maximum theta
504 if (ConstPropertyExists("MR_THETAMAX")) RemoveConstProperty("MR_THETAMAX");
505
506 // Adds a new maximum theta
507 AddConstProperty("MR_THETAMAX", theta_max);
508
509 //G4cout << "theta_max: " << theta_max << G4endl;
510
511 // Removes an existing minimum energy
512 if (ConstPropertyExists("MR_EMIN")) RemoveConstProperty("MR_EMIN");
513
514 // Adds a new minimum energy
515 AddConstProperty("MR_EMIN", E_min);
516
517 //G4cout << "Emin: " << E_min << G4endl;
518
519 // Removes an existing maximum energy
520 if (ConstPropertyExists("MR_EMAX")) RemoveConstProperty("MR_EMAX");
521
522 // Adds a new maximum energy
523 AddConstProperty("MR_EMAX", E_max);
524
525 //G4cout << "Emax: " << E_max << G4endl;
526
527 // Removes an existing Theta angle number
528 if(ConstPropertyExists("MR_ANGNOTHETA"))RemoveConstProperty("MR_ANGNOTHETA");
529
530 // Adds a new Theta angle number
531 AddConstProperty("MR_ANGNOTHETA", (G4double)AngNoTheta);
532
533 //G4cout << "AngNoTheta: " << AngNoTheta << G4endl;
534
535 // Removes an existing Phi angle number
536 if (ConstPropertyExists("MR_ANGNOPHI")) RemoveConstProperty("MR_ANGNOPHI");
537
538 // Adds a new Phi angle number
539 AddConstProperty("MR_ANGNOPHI", (G4double)AngNoPhi);
540
541 //G4cout << "AngNoPhi: " << AngNoPhi << G4endl;
542
543 // Removes an existing angular cut
544 if (ConstPropertyExists("MR_ANGCUT")) RemoveConstProperty("MR_ANGCUT");
545
546 // Adds a new angle number
547 AddConstProperty("MR_ANGCUT", AngularCut);
548
549 //G4cout << "AngularCut: " << AngularCut/degree << "degree" << G4endl;
550
551 // Starts the lookup table calculation
552
554
555 //G4cout << " *** DONE! ***" << G4endl;
556}

References G4MaterialPropertiesTable::AddConstProperty(), ComputeMicroRoughnessTables(), G4MaterialPropertiesTable::ConstPropertyExists(), and G4MaterialPropertiesTable::RemoveConstProperty().

◆ SetMRMaxProbability()

void G4UCNMaterialPropertiesTable::SetMRMaxProbability ( G4double  theta_i,
G4double  Energy,
G4double  value 
)

Definition at line 326 of file G4UCNMaterialPropertiesTable.cc.

328{
330
331 if (theta_i<theta_i_min || theta_i>theta_i_max ||
332 Energy<Emin || Energy>Emax) {
333 } else {
334
335 // Determines the nearest cell in the lookup table which contains
336 // the probability
337
338 G4int theta_i_pos = G4int((theta_i-theta_i_min)/theta_i_step+0.5);
339 G4int E_pos = G4int((Energy-Emin)/E_step+0.5);
340
341 // lookup table is onedimensional (1 row), energy is in rows,
342 // theta_i in columns
343
344 *(maxMicroRoughnessTable+E_pos+theta_i_pos*noE) = value;
345 }
346 }
347}

References E_step, Emax, Emin, maxMicroRoughnessTable, noE, theta_i_max, theta_i_min, and theta_i_step.

◆ SetMRMaxTransProbability()

void G4UCNMaterialPropertiesTable::SetMRMaxTransProbability ( G4double  theta_i,
G4double  Energy,
G4double  value 
)

Definition at line 372 of file G4UCNMaterialPropertiesTable.cc.

374{
376
377 if (theta_i<theta_i_min || theta_i>theta_i_max ||
378 Energy<Emin || Energy>Emax) {
379 } else {
380
381 // Determines the nearest cell in the lookup table which contains
382 // the probability
383
384 G4int theta_i_pos = G4int((theta_i-theta_i_min)/theta_i_step+0.5);
385 G4int E_pos = G4int((Energy-Emin)/E_step+0.5);
386
387 // lookup table is onedimensional (1 row), energy is in rows,
388 // theta_i in columns
389
390 *(maxMicroRoughnessTransTable+E_pos+theta_i_pos*noE) = value;
391 }
392 }
393}

References E_step, Emax, Emin, maxMicroRoughnessTransTable, noE, theta_i_max, theta_i_min, and theta_i_step.

◆ TransConditionsValid()

G4bool G4UCNMaterialPropertiesTable::TransConditionsValid ( G4double  E,
G4double  VFermi,
G4double  theta_i 
)

Definition at line 433 of file G4UCNMaterialPropertiesTable.cc.

436{
438 G4double k_l2 = 2*neutron_mass_c2*VFermi / hbarc_squared;
439
440 if (E*(std::cos(theta_i)*std::cos(theta_i)) < VFermi) return false;
441
442 G4double kS2 = k_l2 - k2;
443
444 //G4cout << "Conditions; 2bk' cos(theta): " << 2*b*sqrt(kS2)*cos(theta_i) <<
445 // ", 2bk_l: " << 2*b*sqrt(k_l2) << G4endl;
446
447 // see eq. 18 of the Steyerl paper
448
449 if (2*b*std::sqrt(kS2)*std::cos(theta_i) < 1 && 2*b*std::sqrt(k_l2) < 1) return true;
450 else return false;
451}

References b, source.hepunit::hbarc_squared, and source.hepunit::neutron_mass_c2.

Field Documentation

◆ AngCut

G4double G4UCNMaterialPropertiesTable::AngCut
private

◆ b

G4double G4UCNMaterialPropertiesTable::b
private

◆ E_step

G4double G4UCNMaterialPropertiesTable::E_step
private

◆ Emax

G4double G4UCNMaterialPropertiesTable::Emax
private

◆ Emin

G4double G4UCNMaterialPropertiesTable::Emin
private

◆ fMatConstPropNames

std::vector<G4String> G4MaterialPropertiesTable::fMatConstPropNames
privateinherited

◆ fMatPropNames

std::vector<G4String> G4MaterialPropertiesTable::fMatPropNames
privateinherited

◆ fMCP

std::vector<std::pair<G4double, G4bool> > G4MaterialPropertiesTable::fMCP
privateinherited

◆ fMP

std::vector<G4MaterialPropertyVector*> G4MaterialPropertiesTable::fMP
privateinherited

◆ maxMicroRoughnessTable

G4double* G4UCNMaterialPropertiesTable::maxMicroRoughnessTable
private

◆ maxMicroRoughnessTransTable

G4double* G4UCNMaterialPropertiesTable::maxMicroRoughnessTransTable
private

◆ no_theta_i

G4int G4UCNMaterialPropertiesTable::no_theta_i
private

◆ noE

G4int G4UCNMaterialPropertiesTable::noE
private

◆ theMicroRoughnessTable

G4double* G4UCNMaterialPropertiesTable::theMicroRoughnessTable
private

◆ theMicroRoughnessTransTable

G4double* G4UCNMaterialPropertiesTable::theMicroRoughnessTransTable
private

◆ theta_i_max

G4double G4UCNMaterialPropertiesTable::theta_i_max
private

◆ theta_i_min

G4double G4UCNMaterialPropertiesTable::theta_i_min
private

◆ theta_i_step

G4double G4UCNMaterialPropertiesTable::theta_i_step
private

◆ w

G4double G4UCNMaterialPropertiesTable::w
private

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