Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PenelopeRayleighModelMI.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id: G4PenelopeRayleighModelMI.hh 75573 2013-11-04 11:48:15Z gcosmo $
27//
28// Author: Luciano Pandola and Gianfranco Paternò
29//
30// -------------------------------------------------------------------
31// History:
32// 03 Dec 2009 L. Pandola 1st implementation
33// 25 May 2011 L. Pandola Renamed (make v2008 as default Penelope)
34// 27 Sep 2013 L. Pandola Migration to MT paradigm
35// 20 Aug 2017 G. Paternò Molecular Interference implementation
36// 24 Mar 2019 G. Paternò Improved Molecular Interference implementation
37// 20 Jun 2020 G. Paternò Read qext separately and leave original atomic form factors
38// 27 Aug 2020 G. Paternò Further improvement of MI implementation
39//
40// -------------------------------------------------------------------
41// Class description:
42// Low Energy Electromagnetic Physics, Rayleigh Scattering
43// with the model from Penelope, version 2008
44// -------------------------------------------------------------------
45//
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47
49
54#include "G4DynamicParticle.hh"
55#include "G4ElementTable.hh"
56#include "G4Element.hh"
58#include "G4AutoLock.hh"
59#include "G4Exp.hh"
60#include "G4ExtendedMaterial.hh"
61#include "G4CrystalExtension.hh"
62#include "G4EmParameters.hh"
63
65#include "G4SystemOfUnits.hh"
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67
69 const G4String& nam) :
70 G4VEmModel(nam),
71 fParticleChange(nullptr),fParticle(nullptr),isInitialised(false),logAtomicCrossSection(nullptr),
72 atomicFormFactor(nullptr),MolInterferenceData(nullptr),logFormFactorTable(nullptr),
73 pMaxTable(nullptr),samplingTable(nullptr),fLocalTable(false),fIsMIActive(true),
74 angularFunction(nullptr), fKnownMaterials(nullptr)
75{
76 fIntrinsicLowEnergyLimit = 100.0*eV;
77 fIntrinsicHighEnergyLimit = 100.0*GeV;
78 //SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
79 SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
80
81 if (part) SetParticle(part);
82
83 verboseLevel = 0;
84 // Verbosity scale:
85 // 0 = nothing
86 // 1 = warning for energy non-conservation
87 // 2 = details of energy budget
88 // 3 = calculation of FF and CS, file openings, sampling of atoms
89 // 4 = entering in methods
90
91 //build the energy grid. It is the same for all materials
92 G4double logenergy = G4Log(fIntrinsicLowEnergyLimit/2.);
93 G4double logmaxenergy = G4Log(1.5*fIntrinsicHighEnergyLimit);
94 //finer grid below 160 keV
95 G4double logtransitionenergy = G4Log(160*keV);
96 G4double logfactor1 = G4Log(10.)/250.;
97 G4double logfactor2 = logfactor1*10;
98 logEnergyGridPMax.push_back(logenergy);
99 do {
100 if (logenergy < logtransitionenergy)
101 logenergy += logfactor1;
102 else
103 logenergy += logfactor2;
104 logEnergyGridPMax.push_back(logenergy);
105 } while (logenergy < logmaxenergy);
106}
107
108//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
109
111{
112 if (IsMaster() || fLocalTable) {
113
114 if (logAtomicCrossSection) {
115 for (auto& item : (*logAtomicCrossSection))
116 if (item.second) delete item.second;
117 delete logAtomicCrossSection;
118 logAtomicCrossSection = nullptr;
119 }
120
121 if (atomicFormFactor) {
122 for (auto& item : (*atomicFormFactor))
123 if (item.second) delete item.second;
124 delete atomicFormFactor;
125 atomicFormFactor = nullptr;
126 }
127
128 if (MolInterferenceData) {
129 for (auto& item : (*MolInterferenceData))
130 if (item.second) delete item.second;
131 delete MolInterferenceData;
132 MolInterferenceData = nullptr;
133 }
134
135 if (fKnownMaterials)
136 {
137 delete fKnownMaterials;
138 fKnownMaterials = nullptr;
139 }
140
141 if (angularFunction)
142 {
143 delete angularFunction;
144 angularFunction = nullptr;
145 }
146
147 ClearTables();
148
149 }
150}
151
152//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
153
154void G4PenelopeRayleighModelMI::ClearTables()
155{
156 /*
157 if (!IsMaster())
158 //Should not be here!
159 G4Exception("G4PenelopeRayleighModelMI::ClearTables()",
160 "em0100",FatalException,"Worker thread in this method");
161 */
162
163 if (logFormFactorTable) {
164 for (auto& item : (*logFormFactorTable))
165 if (item.second) delete item.second;
166 delete logFormFactorTable;
167 logFormFactorTable = nullptr; //zero explicitly
168 }
169
170 if (pMaxTable) {
171 for (auto& item : (*pMaxTable))
172 if (item.second) delete item.second;
173 delete pMaxTable;
174 pMaxTable = nullptr; //zero explicitly
175 }
176
177 if (samplingTable) {
178 for (auto& item : (*samplingTable))
179 if (item.second) delete item.second;
180 delete samplingTable;
181 samplingTable = nullptr; //zero explicitly
182 }
183
184 return;
185}
186
187//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
188
190 const G4DataVector& )
191{
192 if (verboseLevel > 3)
193 G4cout << "Calling G4PenelopeRayleighModelMI::Initialise()" << G4endl;
194
195 SetParticle(part);
196
197
198 if (verboseLevel)
199 G4cout << "# Molecular Interference is " << (fIsMIActive ? "ON" : "OFF") << " #" << G4endl;
200
201 //Only the master model creates/fills/destroys the tables
202 if (IsMaster() && part == fParticle) {
203 //clear tables depending on materials, not the atomic ones
204 ClearTables();
205
206 //Use here the highest verbosity, from G4EmParameter or local
207 G4int globVerb = G4EmParameters::Instance()->Verbose();
208 if (globVerb > verboseLevel)
209 {
210 verboseLevel = globVerb;
211 if (verboseLevel)
212 G4cout << "Verbosity level of G4PenelopeRayleighModelMI set to " << verboseLevel <<
213 " from G4EmParameters()" << G4endl;
214 }
215
216 if (verboseLevel > 3)
217 G4cout << "Calling G4PenelopeRayleighModelMI::Initialise() [master]" << G4endl;
218
219 //Load the list of known materials and the DCS integration grid
220 if (fIsMIActive)
221 {
222 if (!fKnownMaterials)
223 fKnownMaterials = new std::map<G4String,G4String>;
224 if (!fKnownMaterials->size())
225 LoadKnownMIFFMaterials();
226 if (!angularFunction)
227 {
228 //Create and fill once
229 angularFunction = new G4PhysicsFreeVector(Ntheta);
230 CalculateThetaAndAngFun(); //angular funtion for DCS integration
231 }
232 }
233
234 //create new tables
235 //logAtomicCrossSection and atomicFormFactor are created only once,
236 //since they are never cleared
237 if (!logAtomicCrossSection)
238 logAtomicCrossSection = new std::map<G4int,G4PhysicsFreeVector*>;
239
240 if (!atomicFormFactor)
241 atomicFormFactor = new std::map<G4int,G4PhysicsFreeVector*>;
242
243 if (fIsMIActive && !MolInterferenceData)
244 MolInterferenceData = new std::map<G4String,G4PhysicsFreeVector*>; //G. Paternò
245
246 if (!logFormFactorTable)
247 logFormFactorTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
248
249 if (!pMaxTable)
250 pMaxTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
251
252 if (!samplingTable)
253 samplingTable = new std::map<const G4Material*,G4PenelopeSamplingData*>;
254
255 //loop on the used materials
257
258 for (size_t i=0;i<theCoupleTable->GetTableSize();i++) {
259 const G4Material* material =
260 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
261 const G4ElementVector* theElementVector = material->GetElementVector();
262
263 for (size_t j=0;j<material->GetNumberOfElements();j++) {
264 G4int iZ = theElementVector->at(j)->GetZasInt();
265 //read data files only in the master
266 if (!logAtomicCrossSection->count(iZ))
267 ReadDataFile(iZ);
268 }
269
270 //1) Read MI form factors
271 if (fIsMIActive && !MolInterferenceData->count(material->GetName()))
272 ReadMolInterferenceData(material->GetName());
273
274 //2) If the table has not been built for the material, do it!
275 if (!logFormFactorTable->count(material))
276 BuildFormFactorTable(material);
277
278 //3) retrieve or build the sampling table
279 if (!(samplingTable->count(material)))
280 InitializeSamplingAlgorithm(material);
281
282 //4) retrieve or build the pMax data
283 if (!pMaxTable->count(material))
284 GetPMaxTable(material);
285 }
286
287 if (verboseLevel > 1) {
288 G4cout << G4endl << "Penelope Rayleigh model v2008 is initialized" << G4endl
289 << "Energy range: "
290 << LowEnergyLimit() / keV << " keV - "
291 << HighEnergyLimit() / GeV << " GeV"
292 << G4endl;
293 }
294
295 }
296
297 if (isInitialised)
298 return;
299 fParticleChange = GetParticleChangeForGamma();
300 isInitialised = true;
301}
302
303//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
304
306 G4VEmModel *masterModel)
307{
308 if (verboseLevel > 3)
309 G4cout << "Calling G4PenelopeRayleighModelMI::InitialiseLocal()" << G4endl;
310
311 //Check that particle matches: one might have multiple master models
312 //(e.g. for e+ and e-)
313 if (part == fParticle) {
314
315 //Get the const table pointers from the master to the workers
316 const G4PenelopeRayleighModelMI* theModel =
317 static_cast<G4PenelopeRayleighModelMI*> (masterModel);
318
319 //Copy pointers to the data tables
320 logAtomicCrossSection = theModel->logAtomicCrossSection;
321 atomicFormFactor = theModel->atomicFormFactor;
322 MolInterferenceData = theModel->MolInterferenceData;
323 logFormFactorTable = theModel->logFormFactorTable;
324 pMaxTable = theModel->pMaxTable;
325 samplingTable = theModel->samplingTable;
326 fKnownMaterials = theModel->fKnownMaterials;
327 angularFunction = theModel->angularFunction;
328
329 //Copy the G4DataVector with the grid
330 logQSquareGrid = theModel->logQSquareGrid;
331
332 //Same verbosity for all workers, as the master
333 verboseLevel = theModel->verboseLevel;
334
335 }
336
337 return;
338}
339
340//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
341
342namespace {G4Mutex PenelopeRayleighModelMutex = G4MUTEX_INITIALIZER;}
344 G4double energy,
345 G4double Z,
346 G4double,
347 G4double,
348 G4double)
349{
350 //Cross section of Rayleigh scattering in Penelope v2008 is calculated by the EPDL97
351 //tabulation, Cuellen et al. (1997), with non-relativistic form factors from Hubbel
352 //et al. J. Phys. Chem. Ref. Data 4 (1975) 471; Erratum ibid. 6 (1977) 615.
353
354 if (verboseLevel > 3)
355 G4cout << "Calling CrossSectionPerAtom() of G4PenelopeRayleighModelMI" << G4endl;
356
357 G4int iZ = (G4int) Z;
358
359 //Either Initialize() was not called, or we are in a slave and InitializeLocal() was
360 //not invoked
361 if (!logAtomicCrossSection) {
362 //create a **thread-local** version of the table. Used only for G4EmCalculator and
363 //Unit Tests
364 fLocalTable = true;
365 logAtomicCrossSection = new std::map<G4int,G4PhysicsFreeVector*>;
366 }
367
368 //now it should be ok
369 if (!logAtomicCrossSection->count(iZ)) {
370 //If we are here, it means that Initialize() was inkoved, but the MaterialTable was
371 //not filled up. This can happen in a UnitTest or via G4EmCalculator
372 if (verboseLevel > 0) {
373 //Issue a G4Exception (warning) only in verbose mode
375 ed << "Unable to retrieve the cross section table for Z=" << iZ << G4endl;
376 ed << "This can happen only in Unit Tests or via G4EmCalculator" << G4endl;
377 G4Exception("G4PenelopeRayleighModelMI::ComputeCrossSectionPerAtom()",
378 "em2040",JustWarning,ed);
379 }
380
381 //protect file reading via autolock
382 G4AutoLock lock(&PenelopeRayleighModelMutex);
383 ReadDataFile(iZ);
384 lock.unlock();
385
386 }
387
388 G4double cross = 0;
389
390 G4PhysicsFreeVector* atom = logAtomicCrossSection->find(iZ)->second;
391 if (!atom) {
393 ed << "Unable to find Z=" << iZ << " in the atomic cross section table" << G4endl;
394 G4Exception("G4PenelopeRayleighModelMI::ComputeCrossSectionPerAtom()",
395 "em2041",FatalException,ed);
396 return 0;
397 }
398
399 G4double logene = G4Log(energy);
400 G4double logXS = atom->Value(logene);
401 cross = G4Exp(logXS);
402
403 if (verboseLevel > 2) {
404 G4cout << "Rayleigh cross section at " << energy/keV << " keV for Z="
405 << Z << " = " << cross/barn << " barn" << G4endl;
406 }
407
408 return cross;
409}
410
411//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
412
413void G4PenelopeRayleighModelMI::CalculateThetaAndAngFun()
414{
415 G4double theta = 0;
416
417 for(G4int k=0; k<Ntheta; k++) {
418 theta += fDTheta;
419 G4double value = (1+std::cos(theta)*std::cos(theta))*std::sin(theta);
420 angularFunction->PutValue(k,theta,value);
421 if (verboseLevel > 3)
422 G4cout << "theta[" << k << "]: " << angularFunction->Energy(k)
423 << ", angFun[" << k << "]: " << (*angularFunction)[k] << G4endl;
424 }
425}
426
427//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
428
429G4double G4PenelopeRayleighModelMI::CalculateQSquared(G4double angle, G4double energy)
430{
431 G4double lambda,x,q,q2 = 0;
432
433 lambda = hbarc*twopi/energy;
434 x = 1./lambda*std::sin(angle/2.);
435 q = 2.*h_Planck*x/(electron_mass_c2/c_light);
436
437 if (verboseLevel > 3) {
438 G4cout << "E: " << energy/keV << " keV, lambda: " << lambda/nm << " nm"
439 << ", x: " << x*nm << ", q: " << q << G4endl;
440 }
441
442 q2 = std::pow(q,2);
443
444 return q2;
445}
446
447//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
448
449//Overriding of parent's (G4VEmModel) method
451 const G4ParticleDefinition* p,
452 G4double energy,
453 G4double,
454 G4double)
455{
456 //check if we are in a Unit Test (only for the first time)
457 static G4bool amInAUnitTest = false;
458 if (G4ProductionCutsTable::GetProductionCutsTable()->GetTableSize() == 0 && !amInAUnitTest)
459 {
460 amInAUnitTest = true;
462 ed << "The ProductionCuts table is empty " << G4endl;
463 ed << "This should happen only in Unit Tests" << G4endl;
464 G4Exception("G4PenelopeRayleighModelMI::CrossSectionPerVolume()",
465 "em2019",JustWarning,ed);
466 }
467 //If the material does not have a MIFF, continue with the old-style calculation
468 G4String matname = material->GetName();
469 if (amInAUnitTest)
470 {
471 //No need to lock, as this is always called in a master
472 const G4ElementVector* theElementVector = material->GetElementVector();
473 //protect file reading via autolock
474 for (size_t j=0;j<material->GetNumberOfElements();j++) {
475 G4int iZ = theElementVector->at(j)->GetZasInt();
476 if (!logAtomicCrossSection->count(iZ)) {
477 ReadDataFile(iZ);
478 }
479 }
480 if (fIsMIActive)
481 ReadMolInterferenceData(matname);
482 if (!logFormFactorTable->count(material))
483 BuildFormFactorTable(material);
484 if (!(samplingTable->count(material)))
485 InitializeSamplingAlgorithm(material);
486 if (!pMaxTable->count(material))
487 GetPMaxTable(material);
488 }
489
490 G4bool useMIFF = fIsMIActive && (MolInterferenceData->count(matname) || matname.find("MedMat") != std::string::npos);
491 if (!useMIFF)
492 {
493 if (verboseLevel > 2)
494 G4cout << "Rayleigh CS of: " << matname << " calculated through CSperAtom!" << G4endl;
495 return G4VEmModel::CrossSectionPerVolume(material,p,energy);
496 }
497
498 // If we are here, it means that we have to integrate the cross section
499 if (verboseLevel > 2)
500 G4cout << "Rayleigh CS of: " << matname
501 << " calculated through integration of the DCS" << G4endl;
502
503
504 G4double cs = 0;
505
506 //force null cross-section if below the low-energy edge of the table
507 if (energy < LowEnergyLimit()) return cs;
508
509
510
511
512 //if the material is a CRYSTAL, forbid this process
513 if (material->IsExtended() && material->GetName() != "CustomMat") {
514 G4ExtendedMaterial* extendedmaterial = (G4ExtendedMaterial*)material;
515 G4CrystalExtension* crystalExtension = (G4CrystalExtension*)extendedmaterial->RetrieveExtension("crystal");
516 if (crystalExtension != 0) {
517 G4cout << "The material has a crystalline structure, a dedicated diffraction model is used!" << G4endl;
518 return 0;
519 }
520 }
521
522
523 //GET MATERIAL INFORMATION
524
525 G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
526 G4int nElements = material->GetNumberOfElements();
527 const G4ElementVector* elementVector = material->GetElementVector();
528 const G4double* fractionVector = material->GetFractionVector();
529 //const G4double* theAtomNumDensityVector = material-> GetVecNbOfAtomsPerVolume();
530
531 //Stoichiometric factors
532 std::vector<G4double> *StoichiometricFactors = new std::vector<G4double>;
533 for (G4int i=0;i<nElements;i++) {
534 G4double fraction = fractionVector[i];
535 G4double atomicWeigth = (*elementVector)[i]->GetA()/(g/mole);
536 StoichiometricFactors->push_back(fraction/atomicWeigth);
537 }
538 G4double MaxStoichiometricFactor = 0.;
539 for (G4int i=0;i<nElements;i++) {
540 if ((*StoichiometricFactors)[i] > MaxStoichiometricFactor)
541 MaxStoichiometricFactor = (*StoichiometricFactors)[i];
542 }
543 for (G4int i=0;i<nElements;i++) {
544 (*StoichiometricFactors)[i] /= MaxStoichiometricFactor;
545 //G4cout << "StoichiometricFactors[" << i+1 << "]: " << (*StoichiometricFactors)[i] << G4endl;
546 }
547
548 //Equivalent atoms per molecule
549 G4double atPerMol = 0;
550 for (G4int i=0;i<nElements;i++)
551 atPerMol += (*StoichiometricFactors)[i];
552 G4double moleculeDensity = 0.;
553 if (atPerMol) moleculeDensity = atomDensity/atPerMol;
554
555 if (verboseLevel > 2)
556 G4cout << "Material " << material->GetName() << " has " << atPerMol << " atoms "
557 << "per molecule and " << moleculeDensity/(cm*cm*cm) << " molecule/cm3" << G4endl;
558
559 //Equivalent molecular weight (dimensionless)
560 G4double MolWeight = 0.;
561 for (G4int i=0;i<nElements;i++)
562 MolWeight += (*StoichiometricFactors)[i]*(*elementVector)[i]->GetA()/(g/mole);
563
564 if (verboseLevel > 2)
565 G4cout << "Molecular weight of " << matname << ": " << MolWeight << " g/mol" << G4endl;
566
567 G4double IntegrandFun[Ntheta];
568 for (G4int k=0; k<Ntheta; k++) {
569 G4double theta = angularFunction->Energy(k); //the x-value is called "Energy", but is an angle
570 G4double F2 = GetFSquared(material,CalculateQSquared(theta,energy));
571 IntegrandFun[k] = (*angularFunction)[k]*F2;
572 }
573
574 G4double constant = pi*classic_electr_radius*classic_electr_radius;
575 cs = constant*IntegrateFun(IntegrandFun,Ntheta,fDTheta);
576
577 //Now cs is the cross section per molecule, let's calculate the cross section per volume
578 G4double csvolume = cs*moleculeDensity;
579
580
581 //print CS and mfp
582 if (verboseLevel > 2)
583 G4cout << "Rayleigh CS of " << matname << " at " << energy/keV
584 << " keV: " << cs/barn << " barn"
585 << ", mean free path: " << 1./csvolume/mm << " mm" << G4endl;
586
587
588 delete StoichiometricFactors;
589
590 //return CS
591 return csvolume;
592}
593
594//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
595
596void G4PenelopeRayleighModelMI::BuildFormFactorTable(const G4Material* material)
597{
598 if (verboseLevel > 3)
599 G4cout << "Calling BuildFormFactorTable() of G4PenelopeRayleighModelMI" << G4endl;
600
601
602 //GET MATERIAL INFORMATION
603 G4int nElements = material->GetNumberOfElements();
604 const G4ElementVector* elementVector = material->GetElementVector();
605 const G4double* fractionVector = material->GetFractionVector();
606
607 //Stoichiometric factors
608 std::vector<G4double> *StoichiometricFactors = new std::vector<G4double>;
609 for (G4int i=0;i<nElements;i++) {
610 G4double fraction = fractionVector[i];
611 G4double atomicWeigth = (*elementVector)[i]->GetA()/(g/mole);
612 StoichiometricFactors->push_back(fraction/atomicWeigth);
613 }
614 G4double MaxStoichiometricFactor = 0.;
615 for (G4int i=0;i<nElements;i++) {
616 if ((*StoichiometricFactors)[i] > MaxStoichiometricFactor)
617 MaxStoichiometricFactor = (*StoichiometricFactors)[i];
618 }
619 if (MaxStoichiometricFactor<1e-16) {
621 ed << "Inconsistent data of atomic composition for " << material->GetName() << G4endl;
622 G4Exception("G4PenelopeRayleighModelMI::BuildFormFactorTable()",
623 "em2042",FatalException,ed);
624 }
625 for (G4int i=0;i<nElements;i++)
626 (*StoichiometricFactors)[i] /= MaxStoichiometricFactor;
627
628 //Equivalent atoms per molecule
629 G4double atomsPerMolecule = 0;
630 for (G4int i=0;i<nElements;i++)
631 atomsPerMolecule += (*StoichiometricFactors)[i];
632
633 //Equivalent molecular weight (dimensionless)
634 G4double MolWeight = 0.;
635 for (G4int i=0;i<nElements;i++)
636 MolWeight += (*StoichiometricFactors)[i]*(*elementVector)[i]->GetA()/(g/mole);
637
638
639 //CREATE THE FORM FACTOR TABLE
640 //First, the form factors are retrieved [F/sqrt(W)].
641 //Then, they are squared and multiplied for MolWeight -> F2 [dimensionless].
642 //This makes difference for CS calculation, but not for theta sampling.
643 G4PhysicsFreeVector* theFFVec = new G4PhysicsFreeVector(logQSquareGrid.size());
644 theFFVec->SetSpline(true);
645
646 G4String matname = material->GetName();
647 G4String aFileNameFF = "";
648
649 //retrieve MIdata (fFileNameFF)
650 G4MIData* dataMI = GetMIData(material);
651 if (dataMI)
652 aFileNameFF = dataMI->GetFilenameFF();
653
654 //read the MIFF from a file passed by the user
655 if (fIsMIActive && aFileNameFF != "") {
656 if (verboseLevel)
657 G4cout << "Read MIFF for " << matname << " from custom file: " << aFileNameFF << G4endl;
658
659 ReadMolInterferenceData(matname,aFileNameFF);
660 G4PhysicsFreeVector* Fvector = MolInterferenceData->find(matname)->second;
661
662 for (size_t k=0;k<logQSquareGrid.size();k++) {
663 G4double ff2 = 0;
664
665 G4double q = std::pow(G4Exp(logQSquareGrid[k]),0.5);
666 G4double f = Fvector->Value(q);
667 //G4cout << "QGrid[" << k+1 << "]: " << q << ", F[" << k+1 << "]: " << f << G4endl;
668
669 ff2 = f*f*MolWeight;
670 if (ff2) theFFVec->PutValue(k,logQSquareGrid[k],G4Log(ff2));
671 }
672
673 }
674 //retrieve the MIFF from the database or use the IAM
675 else {
676 //medical material: composition of fat, water, bonematrix, mineral
677 if (fIsMIActive && matname.find("MedMat") != std::string::npos) {
678 if (verboseLevel)
679 G4cout << "Build MIFF from components for " << matname << G4endl;
680
681 //get the material composition from its name
682 G4int ki, kf=6, ktot=19;
683 G4double comp[4];
684 G4String compstring = matname.substr(kf+1, ktot);
685 for (size_t j=0; j<4; j++) {
686 ki = kf+1;
687 kf = ki+4;
688 compstring = matname.substr(ki, 4);
689 comp[j] = atof(compstring.c_str());
690 if (verboseLevel > 2)
691 G4cout << " -- MedMat comp[" << j+1 << "]: " << comp[j] << G4endl;
692 }
693
694 char* path = std::getenv("G4LEDATA");
695 if (!path) {
696 G4String excep = "G4LEDATA environment variable not set!";
697 G4Exception("G4PenelopeRayleighModelMI::BuildFormFactorTable()",
698 "em0006",FatalException,excep);
699 }
700
701 if (!MolInterferenceData->count("Fat_MI"))
702 ReadMolInterferenceData("Fat_MI");
703 G4PhysicsFreeVector* fatFF = MolInterferenceData->find("Fat_MI")->second;
704
705 if (!MolInterferenceData->count("Water_MI"))
706 ReadMolInterferenceData("Water_MI");
707 G4PhysicsFreeVector* waterFF = MolInterferenceData->find("Water_MI")->second;
708
709 if (!MolInterferenceData->count("BoneMatrix_MI"))
710 ReadMolInterferenceData("BoneMatrix_MI");
711 G4PhysicsFreeVector* bonematrixFF = MolInterferenceData->find("BoneMatrix_MI")->second;
712
713 if (!MolInterferenceData->count("Mineral_MI"))
714 ReadMolInterferenceData("Mineral_MI");
715 G4PhysicsFreeVector* mineralFF = MolInterferenceData->find("Mineral_MI")->second;
716
717 //get and combine the molecular form factors with interference effect
718 for (size_t k=0;k<logQSquareGrid.size();k++) {
719 G4double ff2 = 0;
720 G4double q = std::pow(G4Exp(logQSquareGrid[k]),0.5);
721 G4double ffat = fatFF->Value(q);
722 G4double fwater = waterFF->Value(q);
723 G4double fbonematrix = bonematrixFF->Value(q);
724 G4double fmineral = mineralFF->Value(q);
725 ff2 = comp[0]*ffat*ffat+comp[1]*fwater*fwater+
726 comp[2]*fbonematrix*fbonematrix+comp[3]*fmineral*fmineral;
727 ff2 *= MolWeight;
728 if (ff2) theFFVec->PutValue(k,logQSquareGrid[k],G4Log(ff2));
729 }
730 }
731 //other materials with MIFF (from the database)
732 else if (fIsMIActive && MolInterferenceData->count(matname)) {
733 if (verboseLevel)
734 G4cout << "Read MIFF from database " << matname << G4endl;
735 G4PhysicsFreeVector* FF = MolInterferenceData->find(matname)->second;
736 for (size_t k=0;k<logQSquareGrid.size();k++) {
737 G4double ff2 = 0;
738 G4double q = std::pow(G4Exp(logQSquareGrid[k]),0.5);
739 G4double f = FF->Value(q);
740 ff2 = f*f*MolWeight;
741 if (ff2) theFFVec->PutValue(k,logQSquareGrid[k],G4Log(ff2));
742 }
743
744 }
745
746 //IAM
747 else {
748 if (verboseLevel)
749 G4cout << "FF of " << matname << " calculated according to the IAM" << G4endl;
750 for (size_t k=0;k<logQSquareGrid.size();k++) {
751 G4double ff2 = 0;
752 for (G4int i=0;i<nElements;i++) {
753 G4int iZ = (*elementVector)[i]->GetZasInt();
754 G4PhysicsFreeVector* theAtomVec = atomicFormFactor->find(iZ)->second;
755 G4double q = std::pow(G4Exp(logQSquareGrid[k]),0.5);
756 G4double f = theAtomVec->Value(q);
757 ff2 += f*f*(*StoichiometricFactors)[i];
758 }
759 if (ff2) theFFVec->PutValue(k,logQSquareGrid[k],G4Log(ff2));
760 }
761 }
762 }
763 logFormFactorTable->insert(std::make_pair(material,theFFVec));
764
765 if (verboseLevel > 3)
766 DumpFormFactorTable(material);
767 delete StoichiometricFactors;
768
769 return;
770
771}
772
773//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
774
775void G4PenelopeRayleighModelMI::SampleSecondaries(std::vector<G4DynamicParticle*>*,
776 const G4MaterialCutsCouple* couple,
777 const G4DynamicParticle* aDynamicGamma,
778 G4double,
779 G4double)
780{
781 // Sampling of the Rayleigh final state (namely, scattering angle of the photon)
782 // from the Penelope2008 model. The scattering angle is sampled from the atomic
783 // cross section dOmega/d(cosTheta) from Born ("Atomic Phyisics", 1969), disregarding
784 // anomalous scattering effects. The Form Factor F(Q) function which appears in the
785 // analytical cross section is retrieved via the method GetFSquared(); atomic data
786 // are tabulated for F(Q). Form factor for compounds is calculated according to
787 // the additivity rule. The sampling from the F(Q) is made via a Rational Inverse
788 // Transform with Aliasing (RITA) algorithm; RITA parameters are calculated once
789 // for each material and managed by G4PenelopeSamplingData objects.
790 // The sampling algorithm (rejection method) has efficiency 67% at low energy, and
791 // increases with energy. For E=100 keV the efficiency is 100% and 86% for
792 // hydrogen and uranium, respectively.
793
794 if (verboseLevel > 3)
795 G4cout << "Calling SamplingSecondaries() of G4PenelopeRayleighModelMI" << G4endl;
796
797 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
798
799 if (photonEnergy0 <= fIntrinsicLowEnergyLimit) {
800 fParticleChange->ProposeTrackStatus(fStopAndKill);
801 fParticleChange->SetProposedKineticEnergy(0.);
802 fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
803 return;
804 }
805
806 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
807
808 const G4Material* theMat = couple->GetMaterial();
809
810 //1) Verify if tables are ready
811 //Either Initialize() was not called, or we are in a slave and InitializeLocal() was
812 //not invoked
813 if (!pMaxTable || !samplingTable || !logAtomicCrossSection || !atomicFormFactor ||
814 !logFormFactorTable) {
815 //create a **thread-local** version of the table. Used only for G4EmCalculator and
816 //Unit Tests
817 fLocalTable = true;
818 if (!logAtomicCrossSection)
819 logAtomicCrossSection = new std::map<G4int,G4PhysicsFreeVector*>;
820 if (!atomicFormFactor)
821 atomicFormFactor = new std::map<G4int,G4PhysicsFreeVector*>;
822 if (!logFormFactorTable)
823 logFormFactorTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
824 if (!pMaxTable)
825 pMaxTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
826 if (!samplingTable)
827 samplingTable = new std::map<const G4Material*,G4PenelopeSamplingData*>;
828 if (fIsMIActive && !MolInterferenceData)
829 MolInterferenceData = new std::map<G4String,G4PhysicsFreeVector*>;
830 }
831
832 if (!samplingTable->count(theMat)) {
833 //If we are here, it means that Initialize() was inkoved, but the MaterialTable was
834 //not filled up. This can happen in a UnitTest
835 if (verboseLevel > 0) {
836 //Issue a G4Exception (warning) only in verbose mode
838 ed << "Unable to find the samplingTable data for " <<
839 theMat->GetName() << G4endl;
840 ed << "This can happen only in Unit Tests" << G4endl;
841 G4Exception("G4PenelopeRayleighModelMI::SampleSecondaries()",
842 "em2019",JustWarning,ed);
843 }
844 const G4ElementVector* theElementVector = theMat->GetElementVector();
845 //protect file reading via autolock
846 G4AutoLock lock(&PenelopeRayleighModelMutex);
847 for (size_t j=0;j<theMat->GetNumberOfElements();j++) {
848 G4int iZ = theElementVector->at(j)->GetZasInt();
849 if (!logAtomicCrossSection->count(iZ)) {
850 lock.lock();
851 ReadDataFile(iZ);
852 lock.unlock();
853 }
854 }
855 lock.lock();
856
857 //1) If the table has not been built for the material, do it!
858 if (!logFormFactorTable->count(theMat))
859 BuildFormFactorTable(theMat);
860
861 //2) retrieve or build the sampling table
862 if (!(samplingTable->count(theMat)))
863 InitializeSamplingAlgorithm(theMat);
864
865 //3) retrieve or build the pMax data
866 if (!pMaxTable->count(theMat))
867 GetPMaxTable(theMat);
868
869 lock.unlock();
870 }
871
872 //Ok, restart the job
873 G4PenelopeSamplingData* theDataTable = samplingTable->find(theMat)->second;
874 G4PhysicsFreeVector* thePMax = pMaxTable->find(theMat)->second;
875
876 G4double cosTheta = 1.0;
877
878 //OK, ready to go!
879 G4double qmax = 2.0*photonEnergy0/electron_mass_c2; //this is non-dimensional now
880
881 if (qmax < 1e-10) { //very low momentum transfer
882 G4bool loopAgain=false;
883 do {
884 loopAgain = false;
885 cosTheta = 1.0-2.0*G4UniformRand();
886 G4double G = 0.5*(1+cosTheta*cosTheta);
887 if (G4UniformRand()>G)
888 loopAgain = true;
889 } while(loopAgain);
890 }
891
892 else { //larger momentum transfer
893
894 size_t nData = theDataTable->GetNumberOfStoredPoints();
895 G4double LastQ2inTheTable = theDataTable->GetX(nData-1);
896 G4double q2max = std::min(qmax*qmax,LastQ2inTheTable);
897
898 G4bool loopAgain = false;
899 G4double MaxPValue = thePMax->Value(photonEnergy0);
900 G4double xx=0;
901
902 //Sampling by rejection method. The rejection function is
903 //G = 0.5*(1+cos^2(theta))
904 do {
905 loopAgain = false;
906 G4double RandomMax = G4UniformRand()*MaxPValue;
907 xx = theDataTable->SampleValue(RandomMax);
908
909 //xx is a random value of q^2 in (0,q2max),sampled according to
910 //F(Q^2) via the RITA algorithm
911 if (xx > q2max)
912 loopAgain = true;
913 cosTheta = 1.0-2.0*xx/q2max;
914 G4double G = 0.5*(1+cosTheta*cosTheta);
915 if (G4UniformRand()>G)
916 loopAgain = true;
917 } while(loopAgain);
918 }
919
920 G4double sinTheta = std::sqrt(1-cosTheta*cosTheta);
921
922 //Scattered photon angles. ( Z - axis along the parent photon)
923 G4double phi = twopi * G4UniformRand() ;
924 G4double dirX = sinTheta*std::cos(phi);
925 G4double dirY = sinTheta*std::sin(phi);
926 G4double dirZ = cosTheta;
927
928 //Update G4VParticleChange for the scattered photon
929 G4ThreeVector photonDirection1(dirX, dirY, dirZ);
930
931 photonDirection1.rotateUz(photonDirection0);
932 fParticleChange->ProposeMomentumDirection(photonDirection1) ;
933 fParticleChange->SetProposedKineticEnergy(photonEnergy0) ;
934
935 return;
936}
937
938//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
939
940void G4PenelopeRayleighModelMI::ReadDataFile(const G4int Z)
941{
942
943 if (verboseLevel > 2) {
944 G4cout << "G4PenelopeRayleighModelMI::ReadDataFile()" << G4endl;
945 G4cout << "Going to read Rayleigh data files for Z=" << Z << G4endl;
946 }
947
948 char* path = std::getenv("G4LEDATA");
949 if (!path) {
950 G4String excep = "G4LEDATA environment variable not set!";
951 G4Exception("G4PenelopeRayleighModelMI::ReadDataFile()",
952 "em0006",FatalException,excep);
953 return;
954 }
955
956 //Read first the cross section file (all the files have 250 points)
957 std::ostringstream ost;
958 if (Z>9)
959 ost << path << "/penelope/rayleigh/pdgra" << Z << ".p08";
960 else
961 ost << path << "/penelope/rayleigh/pdgra0" << Z << ".p08";
962 std::ifstream file(ost.str().c_str());
963
964 if (!file.is_open()) {
965 G4String excep = "Data file " + G4String(ost.str()) + " not found!";
966 G4Exception("G4PenelopeRayleighModelMI::ReadDataFile()",
967 "em0003",FatalException,excep);
968 }
969
970 G4int readZ = 0;
971 size_t nPoints = 0;
972 file >> readZ >> nPoints;
973
974 if (readZ != Z || nPoints <= 0 || nPoints >= 5000) {
976 ed << "Corrupted data file for Z=" << Z << G4endl;
977 G4Exception("G4PenelopeRayleighModelMI::ReadDataFile()",
978 "em0005",FatalException,ed);
979 return;
980 }
981
982 G4PhysicsFreeVector* theVec = new G4PhysicsFreeVector((size_t)nPoints);
983 G4double ene=0,f1=0,f2=0,xs=0;
984 for (size_t i=0;i<nPoints;i++) {
985 file >> ene >> f1 >> f2 >> xs;
986 //dimensional quantities
987 ene *= eV;
988 xs *= cm2;
989 theVec->PutValue(i,G4Log(ene),G4Log(xs));
990 if (file.eof() && i != (nPoints-1)) { //file ended too early
992 ed << "Corrupted data file for Z=" << Z << G4endl;
993 ed << "Found less than " << nPoints << " entries" << G4endl;
994 G4Exception("G4PenelopeRayleighModelMI::ReadDataFile()",
995 "em0005",FatalException,ed);
996 }
997 }
998
999 if (!logAtomicCrossSection) {
1000 G4Exception("G4PenelopeRayleighModelMI::ReadDataFile()",
1001 "em2044",FatalException,"Unable to allocate the atomic cross section table");
1002 delete theVec;
1003 return;
1004 }
1005 logAtomicCrossSection->insert(std::make_pair(Z,theVec));
1006 file.close();
1007
1008 //Then, read the extended momentum transfer file
1009 std::ostringstream ost2;
1010 ost2 << path << "/penelope/rayleigh/MIFF/qext.dat";
1011 file.open(ost2.str().c_str());
1012
1013 if (!file.is_open()) {
1014 G4String excep = "Data file " + G4String(ost2.str()) + " not found!";
1015 G4Exception("G4PenelopeRayleighModelMI::ReadDataFile()",
1016 "em0003",FatalException,excep);
1017 }
1018
1019 G4bool fillQGrid = false;
1020 if (!logQSquareGrid.size()) {
1021 fillQGrid = true;
1022 nPoints = 1142;
1023 }
1024
1025 G4double qext = 0;
1026 if (fillQGrid) { //logQSquareGrid filled only one time
1027 for (size_t i=0;i<nPoints;i++) {
1028 file >> qext;
1029 logQSquareGrid.push_back(2.0*G4Log(qext));
1030 }
1031 }
1032 file.close();
1033
1034 //Finally, read the form factor file
1035 std::ostringstream ost3;
1036 if (Z>9)
1037 ost3 << path << "/penelope/rayleigh/pdaff" << Z << ".p08";
1038 else
1039 ost3 << path << "/penelope/rayleigh/pdaff0" << Z << ".p08";
1040 file.open(ost3.str().c_str());
1041
1042 if (!file.is_open()) {
1043 G4String excep = "Data file " + G4String(ost3.str()) + " not found!";
1044 G4Exception("G4PenelopeRayleighModelMI::ReadDataFile()",
1045 "em0003",FatalException,excep);
1046 }
1047
1048 file >> readZ >> nPoints;
1049
1050 if (readZ != Z || nPoints <= 0 || nPoints >= 5000) {
1052 ed << "Corrupted data file for Z=" << Z << G4endl;
1053 G4Exception("G4PenelopeRayleighModelMI::ReadDataFile()",
1054 "em0005",FatalException,ed);
1055 return;
1056 }
1057
1058 G4PhysicsFreeVector* theFFVec = new G4PhysicsFreeVector((size_t)nPoints);
1059 G4double q=0,ff=0,incoh=0;
1060 for (size_t i=0;i<nPoints;i++) {
1061 file >> q >> ff >> incoh; //q and ff are dimensionless (q is in units of (m_e*c))
1062 theFFVec->PutValue(i,q,ff);
1063 if (file.eof() && i != (nPoints-1)) { //file ended too early
1065 ed << "Corrupted data file for Z=" << Z << G4endl;
1066 ed << "Found less than " << nPoints << " entries" << G4endl;
1067 G4Exception("G4PenelopeRayleighModelMI::ReadDataFile()",
1068 "em0005",FatalException,ed);
1069 }
1070 }
1071
1072 if (!atomicFormFactor) {
1073 G4Exception("G4PenelopeRayleighModelMI::ReadDataFile()",
1074 "em2045",FatalException,
1075 "Unable to allocate the atomicFormFactor data table");
1076 delete theFFVec;
1077 return;
1078 }
1079
1080 atomicFormFactor->insert(std::make_pair(Z,theFFVec));
1081 file.close();
1082 return;
1083}
1084
1085//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1086
1087void G4PenelopeRayleighModelMI::ReadMolInterferenceData(const G4String& matname,
1088 const G4String& FFfilename)
1089{
1090
1091 if (verboseLevel > 2) {
1092 G4cout << "G4PenelopeRayleighModelMI::ReadMolInterferenceData() for material " <<
1093 matname << G4endl;
1094 }
1095 G4bool isLocalFile = (FFfilename != "NULL");
1096
1097 char* path = std::getenv("G4LEDATA");
1098 if (!path) {
1099 G4String excep = "G4LEDATA environment variable not set!";
1100 G4Exception("G4PenelopeRayleighModelMI::ReadMolInterferenceData()",
1101 "em0006",FatalException,excep);
1102 return;
1103 }
1104
1105 if (!(fKnownMaterials->count(matname)) && !isLocalFile) //material not found
1106 return;
1107
1108 G4String aFileName = (isLocalFile) ? FFfilename : fKnownMaterials->find(matname)->second;
1109
1110 //if the material has a MIFF, read it from the database
1111 if (aFileName != "") {
1112 if (verboseLevel > 2)
1113 G4cout << "ReadMolInterferenceData(). Read material: " << matname << ", filename: " <<
1114 aFileName << " " <<
1115 (isLocalFile ? "(local)" : "(database)") << G4endl;
1116 std::ifstream file;
1117 std::ostringstream ostIMFF;
1118 if (isLocalFile)
1119 ostIMFF << aFileName;
1120 else
1121 ostIMFF << path << "/penelope/rayleigh/MIFF/" << aFileName;
1122 file.open(ostIMFF.str().c_str());
1123
1124 if (!file.is_open()) {
1125 G4String excep = "Data file " + G4String(ostIMFF.str()) + " not found!";
1126 G4Exception("G4PenelopeRayleighModelMI::ReadMolInterferenceData()",
1127 "em1031",FatalException,excep);
1128 return;
1129 }
1130
1131 //check the number of points in the file
1132 size_t nPoints = 0;
1133 G4double x=0,y=0;
1134 while (file.good()) {
1135 file >> x >> y;
1136 nPoints++;
1137 }
1138 file.close();
1139 nPoints--;
1140 if (verboseLevel > 3)
1141 G4cout << "Number of nPoints: " << nPoints << G4endl;
1142
1143 //read the file
1144 file.open(ostIMFF.str().c_str());
1145 G4PhysicsFreeVector* theFFVec = new G4PhysicsFreeVector((size_t)nPoints);
1146 G4double q=0,ff=0;
1147 for (size_t i=0; i<nPoints; i++) {
1148 file >> q >> ff;
1149
1150 //q and ff are dimensionless (q is in units of (m_e*c))
1151 theFFVec->PutValue(i,q,ff);
1152
1153 //file ended too early
1154 if (file.eof() && i != (nPoints-1)) {
1156 ed << "Corrupted data file" << G4endl;
1157 ed << "Found less than " << nPoints << " entries" << G4endl;
1158 G4Exception("G4PenelopeRayleighModelMI::ReadMolInterferenceData()",
1159 "em1005",FatalException,ed);
1160 }
1161 }
1162 if (!MolInterferenceData) {
1163 G4Exception("G4PenelopeRayleighModelMI::ReadMolInterferenceData()",
1164 "em2145",FatalException,
1165 "Unable to allocate the Molecular Interference data table");
1166 delete theFFVec;
1167 return;
1168 }
1169 file.close();
1170 MolInterferenceData->insert(std::make_pair(matname,theFFVec));
1171 }
1172
1173 return;
1174}
1175
1176//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1177
1178G4double G4PenelopeRayleighModelMI::GetFSquared(const G4Material* mat, const G4double QSquared)
1179{
1180 G4double f2 = 0;
1181 //Input value QSquared could be zero: protect the log() below against
1182 //the FPE exception
1183
1184 //If Q<1e-10, set Q to 1e-10
1185 G4double logQSquared = (QSquared>1e-10) ? G4Log(QSquared) : -23.;
1186 //last value of the table
1187 G4double maxlogQ2 = logQSquareGrid[logQSquareGrid.size()-1];
1188
1189 //now it should be all right
1190 G4PhysicsFreeVector* theVec = logFormFactorTable->find(mat)->second;
1191
1192 if (!theVec) {
1194 ed << "Unable to retrieve F squared table for " << mat->GetName() << G4endl;
1195 G4Exception("G4PenelopeRayleighModelMI::GetFSquared()",
1196 "em2046",FatalException,ed);
1197 return 0;
1198 }
1199
1200 if (logQSquared < -20) { //Q < 1e-9
1201 G4double logf2 = (*theVec)[0]; //first value of the table
1202 f2 = G4Exp(logf2);
1203 }
1204 else if (logQSquared > maxlogQ2)
1205 f2 =0;
1206 else {
1207 //log(Q^2) vs. log(F^2)
1208 G4double logf2 = theVec->Value(logQSquared);
1209 f2 = G4Exp(logf2);
1210 }
1211
1212 if (verboseLevel > 3) {
1213 G4cout << "G4PenelopeRayleighModelMI::GetFSquared() in " << mat->GetName() << G4endl;
1214 G4cout << "Q^2 = " << QSquared << " (units of 1/(m_e*c)); F^2 = " << f2 << G4endl;
1215 }
1216 return f2;
1217}
1218
1219//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1220
1221void G4PenelopeRayleighModelMI::InitializeSamplingAlgorithm(const G4Material* mat)
1222{
1223
1224 G4double q2min = 0;
1225 G4double q2max = 0;
1226 const size_t np = 150; //hard-coded in Penelope
1227 // G4cout << "Init N= " << logQSquareGrid.size() << G4endl;
1228 for (size_t i=1;i<logQSquareGrid.size();i++)
1229 {
1230 G4double Q2 = G4Exp(logQSquareGrid[i]);
1231 if (GetFSquared(mat,Q2) > 1e-35)
1232 {
1233 q2max = G4Exp(logQSquareGrid[i-1]);
1234 }
1235 //G4cout << "Q2= " << Q2 << " q2max= " << q2max << G4endl;
1236 }
1237
1238 size_t nReducedPoints = np/4;
1239
1240 //check for errors
1241 if (np < 16)
1242 {
1243 G4Exception("G4PenelopeRayleighModelMI::InitializeSamplingAlgorithm()",
1244 "em2047",FatalException,
1245 "Too few points to initialize the sampling algorithm");
1246 }
1247 if (q2min > (q2max-1e-10))
1248 {
1249 G4cout << "q2min= " << q2min << " q2max= " << q2max << G4endl;
1250 G4Exception("G4PenelopeRayleighModelMI::InitializeSamplingAlgorithm()",
1251 "em2048",FatalException,
1252 "Too narrow grid to initialize the sampling algorithm");
1253 }
1254
1255 //This is subroutine RITAI0 of Penelope
1256 //Create an object of type G4PenelopeRayleighSamplingData --> store in a map::Material*
1257
1258 //temporary vectors --> Then everything is stored in G4PenelopeSamplingData
1259 G4DataVector* x = new G4DataVector();
1260
1261 /*******************************************************************************
1262 Start with a grid of NUNIF points uniformly spaced in the interval q2min,q2max
1263 ********************************************************************************/
1264 size_t NUNIF = std::min(std::max(((size_t)8),nReducedPoints),np/2);
1265 const G4int nip = 51; //hard-coded in Penelope
1266
1267 G4double dx = (q2max-q2min)/((G4double) NUNIF-1);
1268 x->push_back(q2min);
1269 for (size_t i=1;i<NUNIF-1;i++)
1270 {
1271 G4double app = q2min + i*dx;
1272 x->push_back(app); //increase
1273 }
1274 x->push_back(q2max);
1275
1276 if (verboseLevel> 3)
1277 G4cout << "Vector x has " << x->size() << " points, while NUNIF = " << NUNIF << G4endl;
1278
1279 //Allocate temporary storage vectors
1280 G4DataVector* area = new G4DataVector();
1281 G4DataVector* a = new G4DataVector();
1282 G4DataVector* b = new G4DataVector();
1283 G4DataVector* c = new G4DataVector();
1284 G4DataVector* err = new G4DataVector();
1285
1286 for (size_t i=0;i<NUNIF-1;i++) //build all points but the last
1287 {
1288 //Temporary vectors for this loop
1289 G4DataVector* pdfi = new G4DataVector();
1290 G4DataVector* pdfih = new G4DataVector();
1291 G4DataVector* sumi = new G4DataVector();
1292
1293 G4double dxi = ((*x)[i+1]-(*x)[i])/(G4double (nip-1));
1294 G4double pdfmax = 0;
1295 for (G4int k=0;k<nip;k++)
1296 {
1297 G4double xik = (*x)[i]+k*dxi;
1298 G4double pdfk = std::max(GetFSquared(mat,xik),0.);
1299 pdfi->push_back(pdfk);
1300 pdfmax = std::max(pdfmax,pdfk);
1301 if (k < (nip-1))
1302 {
1303 G4double xih = xik + 0.5*dxi;
1304 G4double pdfIK = std::max(GetFSquared(mat,xih),0.);
1305 pdfih->push_back(pdfIK);
1306 pdfmax = std::max(pdfmax,pdfIK);
1307 }
1308 }
1309
1310 //Simpson's integration
1311 G4double cons = dxi*0.5*(1./3.);
1312 sumi->push_back(0.);
1313 for (G4int k=1;k<nip;k++)
1314 {
1315 G4double previous = (*sumi)[k-1];
1316 G4double next = previous + cons*((*pdfi)[k-1]+4.0*(*pdfih)[k-1]+(*pdfi)[k]);
1317 sumi->push_back(next);
1318 }
1319
1320 G4double lastIntegral = (*sumi)[sumi->size()-1];
1321 area->push_back(lastIntegral);
1322 //Normalize cumulative function
1323 G4double factor = 1.0/lastIntegral;
1324 for (size_t k=0;k<sumi->size();k++)
1325 (*sumi)[k] *= factor;
1326
1327 //When the PDF vanishes at one of the interval end points, its value is modified
1328 if ((*pdfi)[0] < 1e-35)
1329 (*pdfi)[0] = 1e-5*pdfmax;
1330 if ((*pdfi)[pdfi->size()-1] < 1e-35)
1331 (*pdfi)[pdfi->size()-1] = 1e-5*pdfmax;
1332
1333 G4double pli = (*pdfi)[0]*factor;
1334 G4double pui = (*pdfi)[pdfi->size()-1]*factor;
1335 G4double B_temp = 1.0-1.0/(pli*pui*dx*dx);
1336 G4double A_temp = (1.0/(pli*dx))-1.0-B_temp;
1337 G4double C_temp = 1.0+A_temp+B_temp;
1338 if (C_temp < 1e-35)
1339 {
1340 a->push_back(0.);
1341 b->push_back(0.);
1342 c->push_back(1.);
1343 }
1344 else
1345 {
1346 a->push_back(A_temp);
1347 b->push_back(B_temp);
1348 c->push_back(C_temp);
1349 }
1350
1351 //OK, now get ERR(I), the integral of the absolute difference between the rational interpolation
1352 //and the true pdf, extended over the interval (X(I),X(I+1))
1353 G4int icase = 1; //loop code
1354 G4bool reLoop = false;
1355 err->push_back(0.);
1356 do
1357 {
1358 reLoop = false;
1359 (*err)[i] = 0.; //zero variable
1360 for (G4int k=0;k<nip;k++)
1361 {
1362 G4double rr = (*sumi)[k];
1363 G4double pap = (*area)[i]*(1.0+((*a)[i]+(*b)[i]*rr)*rr)*(1.0+((*a)[i]+(*b)[i]*rr)*rr)/
1364 ((1.0-(*b)[i]*rr*rr)*(*c)[i]*((*x)[i+1]-(*x)[i]));
1365 if (k == 0 || k == nip-1)
1366 (*err)[i] += 0.5*std::fabs(pap-(*pdfi)[k]);
1367 else
1368 (*err)[i] += std::fabs(pap-(*pdfi)[k]);
1369 }
1370 (*err)[i] *= dxi;
1371
1372 //If err(I) is too large, the pdf is approximated by a uniform distribution
1373 if ((*err)[i] > 0.1*(*area)[i] && icase == 1)
1374 {
1375 (*b)[i] = 0;
1376 (*a)[i] = 0;
1377 (*c)[i] = 1.;
1378 icase = 2;
1379 reLoop = true;
1380 }
1381 }while(reLoop);
1382
1383 delete pdfi;
1384 delete pdfih;
1385 delete sumi;
1386 } //end of first loop over i
1387
1388 //Now assign last point
1389 (*x)[x->size()-1] = q2max;
1390 a->push_back(0.);
1391 b->push_back(0.);
1392 c->push_back(0.);
1393 err->push_back(0.);
1394 area->push_back(0.);
1395
1396 if (x->size() != NUNIF || a->size() != NUNIF ||
1397 err->size() != NUNIF || area->size() != NUNIF)
1398 {
1400 ed << "Problem in building the Table for Sampling: array dimensions do not match" << G4endl;
1401 G4Exception("G4PenelopeRayleighModelMI::InitializeSamplingAlgorithm()",
1402 "em2049",FatalException,ed);
1403 }
1404
1405 /*******************************************************************************
1406 New grid points are added by halving the sub-intervals with the largest absolute error
1407 This is done up to np=150 points in the grid
1408 ********************************************************************************/
1409 do
1410 {
1411 G4double maxError = 0.0;
1412 size_t iErrMax = 0;
1413 for (size_t i=0;i<err->size()-2;i++)
1414 {
1415 //maxError is the lagest of the interval errors err[i]
1416 if ((*err)[i] > maxError)
1417 {
1418 maxError = (*err)[i];
1419 iErrMax = i;
1420 }
1421 }
1422
1423 //OK, now I have to insert one new point in the position iErrMax
1424 G4double newx = 0.5*((*x)[iErrMax]+(*x)[iErrMax+1]);
1425
1426 x->insert(x->begin()+iErrMax+1,newx);
1427 //Add place-holders in the other vectors
1428 area->insert(area->begin()+iErrMax+1,0.);
1429 a->insert(a->begin()+iErrMax+1,0.);
1430 b->insert(b->begin()+iErrMax+1,0.);
1431 c->insert(c->begin()+iErrMax+1,0.);
1432 err->insert(err->begin()+iErrMax+1,0.);
1433
1434 //Now calculate the other parameters
1435 for (size_t i=iErrMax;i<=iErrMax+1;i++)
1436 {
1437 //Temporary vectors for this loop
1438 G4DataVector* pdfi = new G4DataVector();
1439 G4DataVector* pdfih = new G4DataVector();
1440 G4DataVector* sumi = new G4DataVector();
1441
1442 G4double dxLocal = (*x)[i+1]-(*x)[i];
1443 G4double dxi = ((*x)[i+1]-(*x)[i])/(G4double (nip-1));
1444 G4double pdfmax = 0;
1445 for (G4int k=0;k<nip;k++)
1446 {
1447 G4double xik = (*x)[i]+k*dxi;
1448 G4double pdfk = std::max(GetFSquared(mat,xik),0.);
1449 pdfi->push_back(pdfk);
1450 pdfmax = std::max(pdfmax,pdfk);
1451 if (k < (nip-1))
1452 {
1453 G4double xih = xik + 0.5*dxi;
1454 G4double pdfIK = std::max(GetFSquared(mat,xih),0.);
1455 pdfih->push_back(pdfIK);
1456 pdfmax = std::max(pdfmax,pdfIK);
1457 }
1458 }
1459
1460 //Simpson's integration
1461 G4double cons = dxi*0.5*(1./3.);
1462 sumi->push_back(0.);
1463 for (G4int k=1;k<nip;k++)
1464 {
1465 G4double previous = (*sumi)[k-1];
1466 G4double next = previous + cons*((*pdfi)[k-1]+4.0*(*pdfih)[k-1]+(*pdfi)[k]);
1467 sumi->push_back(next);
1468 }
1469 G4double lastIntegral = (*sumi)[sumi->size()-1];
1470 (*area)[i] = lastIntegral;
1471
1472 //Normalize cumulative function
1473 G4double factor = 1.0/lastIntegral;
1474 for (size_t k=0;k<sumi->size();k++)
1475 (*sumi)[k] *= factor;
1476
1477 //When the PDF vanishes at one of the interval end points, its value is modified
1478 if ((*pdfi)[0] < 1e-35)
1479 (*pdfi)[0] = 1e-5*pdfmax;
1480 if ((*pdfi)[pdfi->size()-1] < 1e-35)
1481 (*pdfi)[pdfi->size()-1] = 1e-5*pdfmax;
1482
1483 G4double pli = (*pdfi)[0]*factor;
1484 G4double pui = (*pdfi)[pdfi->size()-1]*factor;
1485 G4double B_temp = 1.0-1.0/(pli*pui*dxLocal*dxLocal);
1486 G4double A_temp = (1.0/(pli*dxLocal))-1.0-B_temp;
1487 G4double C_temp = 1.0+A_temp+B_temp;
1488 if (C_temp < 1e-35)
1489 {
1490 (*a)[i]= 0.;
1491 (*b)[i] = 0.;
1492 (*c)[i] = 1;
1493 }
1494 else
1495 {
1496 (*a)[i]= A_temp;
1497 (*b)[i] = B_temp;
1498 (*c)[i] = C_temp;
1499 }
1500 //OK, now get ERR(I), the integral of the absolute difference between the rational interpolation
1501 //and the true pdf, extended over the interval (X(I),X(I+1))
1502 G4int icase = 1; //loop code
1503 G4bool reLoop = false;
1504 do
1505 {
1506 reLoop = false;
1507 (*err)[i] = 0.; //zero variable
1508 for (G4int k=0;k<nip;k++)
1509 {
1510 G4double rr = (*sumi)[k];
1511 G4double pap = (*area)[i]*(1.0+((*a)[i]+(*b)[i]*rr)*rr)*(1.0+((*a)[i]+(*b)[i]*rr)*rr)/
1512 ((1.0-(*b)[i]*rr*rr)*(*c)[i]*((*x)[i+1]-(*x)[i]));
1513 if (k == 0 || k == nip-1)
1514 (*err)[i] += 0.5*std::fabs(pap-(*pdfi)[k]);
1515 else
1516 (*err)[i] += std::fabs(pap-(*pdfi)[k]);
1517 }
1518 (*err)[i] *= dxi;
1519
1520 //If err(I) is too large, the pdf is approximated by a uniform distribution
1521 if ((*err)[i] > 0.1*(*area)[i] && icase == 1)
1522 {
1523 (*b)[i] = 0;
1524 (*a)[i] = 0;
1525 (*c)[i] = 1.;
1526 icase = 2;
1527 reLoop = true;
1528 }
1529 }while(reLoop);
1530 delete pdfi;
1531 delete pdfih;
1532 delete sumi;
1533 }
1534 }while(x->size() < np);
1535
1536 if (x->size() != np || a->size() != np ||
1537 err->size() != np || area->size() != np)
1538 {
1539 G4Exception("G4PenelopeRayleighModelMI::InitializeSamplingAlgorithm()",
1540 "em2050",FatalException,
1541 "Problem in building the extended Table for Sampling: array dimensions do not match ");
1542 }
1543
1544 /*******************************************************************************
1545 Renormalization
1546 ********************************************************************************/
1547 G4double ws = 0;
1548 for (size_t i=0;i<np-1;i++)
1549 ws += (*area)[i];
1550 ws = 1.0/ws;
1551 G4double errMax = 0;
1552 for (size_t i=0;i<np-1;i++)
1553 {
1554 (*area)[i] *= ws;
1555 (*err)[i] *= ws;
1556 errMax = std::max(errMax,(*err)[i]);
1557 }
1558
1559 //Vector with the normalized cumulative distribution
1560 G4DataVector* PAC = new G4DataVector();
1561 PAC->push_back(0.);
1562 for (size_t i=0;i<np-1;i++)
1563 {
1564 G4double previous = (*PAC)[i];
1565 PAC->push_back(previous+(*area)[i]);
1566 }
1567 (*PAC)[PAC->size()-1] = 1.;
1568
1569 /*******************************************************************************
1570 Pre-calculated limits for the initial binary search for subsequent sampling
1571 ********************************************************************************/
1572
1573 //G4DataVector* ITTL = new G4DataVector();
1574 std::vector<size_t> *ITTL = new std::vector<size_t>;
1575 std::vector<size_t> *ITTU = new std::vector<size_t>;
1576
1577 //Just create place-holders
1578 for (size_t i=0;i<np;i++)
1579 {
1580 ITTL->push_back(0);
1581 ITTU->push_back(0);
1582 }
1583
1584 G4double bin = 1.0/(np-1);
1585 (*ITTL)[0]=0;
1586 for (size_t i=1;i<(np-1);i++)
1587 {
1588 G4double ptst = i*bin;
1589 G4bool found = false;
1590 for (size_t j=(*ITTL)[i-1];j<np && !found;j++)
1591 {
1592 if ((*PAC)[j] > ptst)
1593 {
1594 (*ITTL)[i] = j-1;
1595 (*ITTU)[i-1] = j;
1596 found = true;
1597 }
1598 }
1599 }
1600 (*ITTU)[ITTU->size()-2] = ITTU->size()-1;
1601 (*ITTU)[ITTU->size()-1] = ITTU->size()-1;
1602 (*ITTL)[ITTL->size()-1] = ITTU->size()-2;
1603
1604 if (ITTU->size() != np || ITTU->size() != np)
1605 {
1606 G4Exception("G4PenelopeRayleighModelMI::InitializeSamplingAlgorithm()",
1607 "em2051",FatalException,
1608 "Problem in building the Limit Tables for Sampling: array dimensions do not match");
1609 }
1610
1611
1612 /********************************************************************************
1613 Copy tables
1614 ********************************************************************************/
1616 for (size_t i=0;i<np;i++)
1617 {
1618 theTable->AddPoint((*x)[i],(*PAC)[i],(*a)[i],(*b)[i],(*ITTL)[i],(*ITTU)[i]);
1619 }
1620
1621 if (verboseLevel > 2)
1622 {
1623 G4cout << "*************************************************************************" <<
1624 G4endl;
1625 G4cout << "Sampling table for Penelope Rayleigh scattering in " << mat->GetName() << G4endl;
1626 theTable->DumpTable();
1627 }
1628 samplingTable->insert(std::make_pair(mat,theTable));
1629
1630
1631 //Clean up temporary vectors
1632 delete x;
1633 delete a;
1634 delete b;
1635 delete c;
1636 delete err;
1637 delete area;
1638 delete PAC;
1639 delete ITTL;
1640 delete ITTU;
1641
1642 //DONE!
1643 return;
1644
1645}
1646
1647//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1648
1649void G4PenelopeRayleighModelMI::GetPMaxTable(const G4Material* mat)
1650{
1651
1652 if (!pMaxTable)
1653 {
1654 G4cout << "G4PenelopeRayleighModelMI::BuildPMaxTable" << G4endl;
1655 G4cout << "Going to instanziate the pMaxTable !" << G4endl;
1656 G4cout << "That should _not_ be here! " << G4endl;
1657 pMaxTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
1658 }
1659 //check if the table is already there
1660 if (pMaxTable->count(mat))
1661 return;
1662
1663 //otherwise build it
1664 if (!samplingTable)
1665 {
1666 G4Exception("G4PenelopeRayleighModelMI::GetPMaxTable()",
1667 "em2052",FatalException,
1668 "SamplingTable is not properly instantiated");
1669 return;
1670 }
1671
1672 //This should not be: the sampling table is built before the p-table
1673 if (!samplingTable->count(mat))
1674 {
1676 ed << "Sampling table for material " << mat->GetName() << " not found";
1677 G4Exception("G4PenelopeRayleighModelMI::GetPMaxTable()",
1678 "em2052",FatalException,
1679 ed);
1680 return;
1681 }
1682
1683 G4PenelopeSamplingData *theTable = samplingTable->find(mat)->second;
1684 size_t tablePoints = theTable->GetNumberOfStoredPoints();
1685
1686 size_t nOfEnergyPoints = logEnergyGridPMax.size();
1687 G4PhysicsFreeVector* theVec = new G4PhysicsFreeVector(nOfEnergyPoints);
1688
1689 const size_t nip = 51; //hard-coded in Penelope
1690
1691 for (size_t ie=0;ie<logEnergyGridPMax.size();ie++)
1692 {
1693 G4double energy = G4Exp(logEnergyGridPMax[ie]);
1694 G4double Qm = 2.0*energy/electron_mass_c2; //this is non-dimensional now
1695 G4double Qm2 = Qm*Qm;
1696 G4double firstQ2 = theTable->GetX(0);
1697 G4double lastQ2 = theTable->GetX(tablePoints-1);
1698 G4double thePMax = 0;
1699
1700 if (Qm2 > firstQ2)
1701 {
1702 if (Qm2 < lastQ2)
1703 {
1704 //bisection to look for the index of Qm
1705 size_t lowerBound = 0;
1706 size_t upperBound = tablePoints-1;
1707 while (lowerBound <= upperBound)
1708 {
1709 size_t midBin = (lowerBound + upperBound)/2;
1710 if( Qm2 < theTable->GetX(midBin))
1711 { upperBound = midBin-1; }
1712 else
1713 { lowerBound = midBin+1; }
1714 }
1715 //upperBound is the output (but also lowerBounf --> should be the same!)
1716 G4double Q1 = theTable->GetX(upperBound);
1717 G4double Q2 = Qm2;
1718 G4double DQ = (Q2-Q1)/((G4double)(nip-1));
1719 G4double theA = theTable->GetA(upperBound);
1720 G4double theB = theTable->GetB(upperBound);
1721 G4double thePAC = theTable->GetPAC(upperBound);
1722 G4DataVector* fun = new G4DataVector();
1723 for (size_t k=0;k<nip;k++)
1724 {
1725 G4double qi = Q1 + k*DQ;
1726 G4double tau = (qi-Q1)/
1727 (theTable->GetX(upperBound+1)-Q1);
1728 G4double con1 = 2.0*theB*tau;
1729 G4double ci = 1.0+theA+theB;
1730 G4double con2 = ci-theA*tau;
1731 G4double etap = 0;
1732 if (std::fabs(con1) > 1.0e-16*std::fabs(con2))
1733 etap = con2*(1.0-std::sqrt(1.0-2.0*tau*con1/(con2*con2)))/con1;
1734 else
1735 etap = tau/con2;
1736 G4double theFun = (theTable->GetPAC(upperBound+1)-thePAC)*
1737 (1.0+(theA+theB*etap)*etap)*(1.0+(theA+theB*etap)*etap)/
1738 ((1.0-theB*etap*etap)*ci*(theTable->GetX(upperBound+1)-Q1));
1739 fun->push_back(theFun);
1740 }
1741 //Now intergrate numerically the fun Cavalieri-Simpson's method
1742 G4DataVector* sum = new G4DataVector;
1743 G4double CONS = DQ*(1./12.);
1744 G4double HCONS = 0.5*CONS;
1745 sum->push_back(0.);
1746 G4double secondPoint = (*sum)[0] +
1747 (5.0*(*fun)[0]+8.0*(*fun)[1]-(*fun)[2])*CONS;
1748 sum->push_back(secondPoint);
1749 for (size_t hh=2;hh<nip-1;hh++)
1750 {
1751 G4double previous = (*sum)[hh-1];
1752 G4double next = previous+(13.0*((*fun)[hh-1]+(*fun)[hh])-
1753 (*fun)[hh+1]-(*fun)[hh-2])*HCONS;
1754 sum->push_back(next);
1755 }
1756 G4double last = (*sum)[nip-2]+(5.0*(*fun)[nip-1]+8.0*(*fun)[nip-2]-
1757 (*fun)[nip-3])*CONS;
1758 sum->push_back(last);
1759 thePMax = thePAC + (*sum)[sum->size()-1]; //last point
1760 delete fun;
1761 delete sum;
1762 }
1763 else
1764 {
1765 thePMax = 1.0;
1766 }
1767 }
1768 else
1769 {
1770 thePMax = theTable->GetPAC(0);
1771 }
1772
1773 //Write number in the table
1774 theVec->PutValue(ie,energy,thePMax);
1775 }
1776
1777 pMaxTable->insert(std::make_pair(mat,theVec));
1778 return;
1779
1780}
1781
1782//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1783
1785{
1786 G4cout << "*****************************************************************" << G4endl;
1787 G4cout << "G4PenelopeRayleighModelMI: Form Factor Table for " << mat->GetName() << G4endl;
1788 //try to use the same format as Penelope-Fortran, namely Q (/m_e*c) and F
1789 G4cout << "Q/(m_e*c) F(Q) " << G4endl;
1790 G4cout << "*****************************************************************" << G4endl;
1791 if (!logFormFactorTable->count(mat))
1792 BuildFormFactorTable(mat);
1793
1794 G4PhysicsFreeVector* theVec = logFormFactorTable->find(mat)->second;
1795 for (size_t i=0;i<theVec->GetVectorLength();i++)
1796 {
1797 G4double logQ2 = theVec->GetLowEdgeEnergy(i);
1798 G4double Q = G4Exp(0.5*logQ2);
1799 G4double logF2 = (*theVec)[i];
1800 G4double F = G4Exp(0.5*logF2);
1801 G4cout << Q << " " << F << G4endl;
1802 }
1803 //DONE
1804 return;
1805}
1806
1807//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...
1808
1809void G4PenelopeRayleighModelMI::SetParticle(const G4ParticleDefinition* p)
1810{
1811 if(!fParticle) {
1812 fParticle = p;
1813 }
1814}
1815
1816//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...
1817void G4PenelopeRayleighModelMI::LoadKnownMIFFMaterials()
1818{
1819 fKnownMaterials->insert(std::pair<G4String,G4String>("Fat_MI","FF_fat_Tartari2002.dat"));
1820 fKnownMaterials->insert(std::pair<G4String,G4String>("Water_MI","FF_water_Tartari2002.dat"));
1821 fKnownMaterials->insert(std::pair<G4String,G4String>("BoneMatrix_MI","FF_bonematrix_Tartari2002.dat"));
1822 fKnownMaterials->insert(std::pair<G4String,G4String>("Mineral_MI","FF_mineral_Tartari2002.dat"));
1823 fKnownMaterials->insert(std::pair<G4String,G4String>("adipose_MI","FF_adipose_Poletti2002.dat"));
1824 fKnownMaterials->insert(std::pair<G4String,G4String>("glandular_MI","FF_glandular_Poletti2002.dat"));
1825 fKnownMaterials->insert(std::pair<G4String,G4String>("breast5050_MI","FF_human_breast_Peplow1998.dat"));
1826 fKnownMaterials->insert(std::pair<G4String,G4String>("carcinoma_MI","FF_carcinoma_Kidane1999.dat"));
1827 fKnownMaterials->insert(std::pair<G4String,G4String>("muscle_MI","FF_pork_muscle_Peplow1998.dat"));
1828 fKnownMaterials->insert(std::pair<G4String,G4String>("kidney_MI","FF_pork_kidney_Peplow1998.dat"));
1829 fKnownMaterials->insert(std::pair<G4String,G4String>("liver_MI","FF_pork_liver_Peplow1998.dat"));
1830 fKnownMaterials->insert(std::pair<G4String,G4String>("heart_MI","FF_pork_heart_Peplow1998.dat"));
1831 fKnownMaterials->insert(std::pair<G4String,G4String>("blood_MI","FF_beef_blood_Peplow1998.dat"));
1832 fKnownMaterials->insert(std::pair<G4String,G4String>("grayMatter_MI","FF_gbrain_DeFelici2008.dat"));
1833 fKnownMaterials->insert(std::pair<G4String,G4String>("whiteMatter_MI","FF_wbrain_DeFelici2008.dat"));
1834 fKnownMaterials->insert(std::pair<G4String,G4String>("bone_MI","FF_bone_King2011.dat"));
1835 fKnownMaterials->insert(std::pair<G4String,G4String>("FatLowX_MI","FF_fat_Tartari2002_joint_lowXdata_ESRF2003.dat"));
1836 fKnownMaterials->insert(std::pair<G4String,G4String>("BoneMatrixLowX_MI","FF_bonematrix_Tartari2002_joint_lowXdata.dat"));
1837 fKnownMaterials->insert(std::pair<G4String,G4String>("PMMALowX_MI", "FF_PMMA_Tartari2002_joint_lowXdata_ESRF2003.dat"));
1838 fKnownMaterials->insert(std::pair<G4String,G4String>("dryBoneLowX_MI","FF_drybone_Tartari2002_joint_lowXdata_ESRF2003.dat"));
1839 fKnownMaterials->insert(std::pair<G4String,G4String>("CIRS30-70_MI","FF_CIRS30-70_Poletti2002.dat"));
1840 fKnownMaterials->insert(std::pair<G4String,G4String>("CIRS50-50_MI","FF_CIRS50-50_Poletti2002.dat"));
1841 fKnownMaterials->insert(std::pair<G4String,G4String>("CIRS70-30_MI","FF_CIRS70-30_Poletti2002.dat"));
1842 fKnownMaterials->insert(std::pair<G4String,G4String>("RMI454_MI", "FF_RMI454_Poletti2002.dat"));
1843 fKnownMaterials->insert(std::pair<G4String,G4String>("PMMA_MI","FF_PMMA_Tartari2002.dat"));
1844 fKnownMaterials->insert(std::pair<G4String,G4String>("Lexan_MI","FF_lexan_Peplow1998.dat"));
1845 fKnownMaterials->insert(std::pair<G4String,G4String>("Kapton_MI","FF_kapton_Peplow1998.dat"));
1846 fKnownMaterials->insert(std::pair<G4String,G4String>("Nylon_MI","FF_nylon_Kosanetzky1987.dat"));
1847 fKnownMaterials->insert(std::pair<G4String,G4String>("Polyethylene_MI","FF_polyethylene_Kosanetzky1987.dat"));
1848 fKnownMaterials->insert(std::pair<G4String,G4String>("Polystyrene_MI","FF_polystyrene_Kosanetzky1987.dat"));
1849 fKnownMaterials->insert(std::pair<G4String,G4String>("Formaline_MI","FF_formaline_Peplow1998.dat"));
1850 fKnownMaterials->insert(std::pair<G4String,G4String>("Acetone_MI","FF_acetone_Cozzini2010.dat"));
1851 fKnownMaterials->insert(std::pair<G4String,G4String>("Hperoxide_MI","FF_Hperoxide_Cozzini2010.dat"));
1852
1853}
1854
1855//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1856
1857G4double G4PenelopeRayleighModelMI::IntegrateFun(G4double y[], G4int n, G4double dTheta)
1858{
1859 G4double integral = 0.;
1860 for (G4int k=0; k<n-1; k++) {
1861 integral += (y[k]+y[k+1]);
1862 }
1863 integral *= dTheta*0.5;
1864 return integral;
1865}
1866
1867//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1868G4MIData* G4PenelopeRayleighModelMI::GetMIData(const G4Material* material)
1869{
1870 if (material->IsExtended()) {
1871 G4ExtendedMaterial* aEM = (G4ExtendedMaterial*)material;
1872 G4MIData* dataMI = (G4MIData*)aEM->RetrieveExtension("MI");
1873 return dataMI;
1874 } else {
1875 return nullptr;
1876 }
1877}
std::vector< G4Element * > G4ElementVector
@ JustWarning
@ 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
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
std::mutex G4Mutex
Definition: G4Threading.hh:81
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4EmParameters * Instance()
G4int Verbose() const
G4VMaterialExtension * RetrieveExtension(const G4String &name)
const G4String & GetFilenameFF()
Definition: G4MIData.hh:53
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:207
virtual G4bool IsExtended() const
Definition: G4Material.cc:795
const G4double * GetFractionVector() const
Definition: G4Material.hh:192
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4String & GetName() const
Definition: G4Material.hh:175
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0., G4double maxEnergy=DBL_MAX) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void DumpFormFactorTable(const G4Material *)
G4PenelopeRayleighModelMI(const G4ParticleDefinition *p=nullptr, const G4String &processName="PenRayleighMI")
G4double GetA(size_t index)
G4double GetPAC(size_t index)
void AddPoint(G4double x0, G4double pac0, G4double a0, G4double b0, size_t ITTL0, size_t ITTU0)
G4double GetX(size_t index)
G4double SampleValue(G4double rndm)
G4double GetB(size_t index)
void PutValue(std::size_t index, G4double energy, G4double dValue)
G4double Energy(std::size_t index) const
G4double Value(G4double theEnergy, std::size_t &lastidx) const
G4double GetLowEdgeEnergy(std::size_t binNumber) const
void SetSpline(G4bool)
std::size_t GetVectorLength() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:254
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double energy(const ThreeVector &p, const G4double m)