Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4NeutronGeneralProcess.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4NeutronGeneralProcess
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 08.08.2022
37//
38// Modifications:
39//
40// Class Description:
41//
42
43// -------------------------------------------------------------------
44//
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47
50#include "G4SystemOfUnits.hh"
51#include "G4HadronicProcess.hh"
53#include "G4Step.hh"
54#include "G4Track.hh"
56#include "G4PhysicsTable.hh"
57#include "G4PhysicsLogVector.hh"
58#include "G4VParticleChange.hh"
61#include "G4Material.hh"
62#include "G4MaterialTable.hh"
63#include "G4Element.hh"
64#include "G4Neutron.hh"
66#include "G4NeutronElasticXS.hh"
67#include "G4NeutronCaptureXS.hh"
68#include "G4Threading.hh"
69
70#include "G4Log.hh"
71#include <iostream>
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74
75G4HadDataHandler* G4NeutronGeneralProcess::theHandler = nullptr;
76
77G4String G4NeutronGeneralProcess::nameT[nTables] = {"0","1","2","3","4"};
78
80: G4HadronicProcess(pname),
81 fMinEnergy(1*CLHEP::keV),
82 fMiddleEnergy(20*CLHEP::MeV),
83 fMaxEnergy(100*CLHEP::TeV),
84 fTimeLimit(10*CLHEP::microsecond)
85{
88
89 fNeutron = G4Neutron::Neutron();
90
92 isMaster = false;
93 }
94}
95
96//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
97
99{
100 if(isMaster) {
101 delete theHandler;
102 theHandler = nullptr;
103 }
104}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107
109{
110 return true;
111}
112
113//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
114
116{
117 fInelastic = ptr;
118 fXSSInelastic = ptr->GetCrossSectionDataStore();
119 fInelasticXS = InitialisationXS(ptr);
120 if(nullptr == fInelasticXS) {
121 fInelasticXS = new G4NeutronInelasticXS();
122 ptr->AddDataSet(fInelasticXS);
123 }
124}
125
126//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127
129{
130 fElastic = ptr;
131 fXSSElastic = ptr->GetCrossSectionDataStore();
132 fElasticXS = InitialisationXS(ptr);
133 if(nullptr == fElasticXS) {
134 fElasticXS = new G4NeutronElasticXS();
135 ptr->AddDataSet(fElasticXS);
136 }
137}
138
139//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
140
142{
143 fCapture = ptr;
144 fXSSCapture = ptr->GetCrossSectionDataStore();
145 fCaptureXS = InitialisationXS(ptr);
146 if(nullptr == fCaptureXS) {
147 fCaptureXS = new G4NeutronCaptureXS();
148 ptr->AddDataSet(fCaptureXS);
149 }
150}
151
152//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
153
155G4NeutronGeneralProcess::InitialisationXS(G4HadronicProcess* proc)
156{
157 G4VCrossSectionDataSet* ptr = nullptr;
158 auto xsv = proc->GetCrossSectionDataStore()->GetDataSetList();
159 if(!xsv.empty()) {
160 ptr = xsv[0];
161 }
162 return ptr;
163}
164
165//....Ooooo0ooooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
166
168{
169 if(1 < verboseLevel) {
170 G4cout << "G4NeutronGeneralProcess::PreparePhysicsTable() for "
171 << GetProcessName()
172 << " and particle " << part.GetParticleName()
173 << " isMaster: " << isMaster << G4endl;
174 }
175 G4bool noEl = (nullptr == fElastic);
176 G4bool noInel = (nullptr == fInelastic);
177 G4bool noCap = (nullptr == fCapture);
178 if(noEl || noInel || noCap) {
180 ed << "Incomplete configuration of the neutron general process." << G4endl;
181 if(noEl) { ed << "Neutron elastic process is not defined" << G4endl; }
182 if(noInel) { ed << "Neutron inelastic process is not defined" << G4endl; }
183 if(noCap) { ed << "Neutron capture process is not defined" << G4endl; }
184 G4Exception ("G4NeutronGeneralProcess::PreparePhysicsTable(..)", "had001",
185 FatalException, ed, "");
186 return;
187 }
188
190
192 fMaxEnergy = std::max(100*MeV, param->GetMaxEnergy());
193 if(param->ApplyFactorXS()) {
194 fXSFactorEl = param->XSFactorNucleonElastic();
195 fXSFactorInel = param->XSFactorNucleonInelastic();
196 }
197
198 fElastic->PreparePhysicsTable(part);
199 fInelastic->PreparePhysicsTable(part);
200 fCapture->PreparePhysicsTable(part);
201
202 std::size_t nmat = G4Material::GetNumberOfMaterials();
204
205 std::size_t nmax = 0;
206 for(std::size_t i=0; i<nmat; ++i) {
207 std::size_t nelm = (*matTable)[i]->GetNumberOfElements();
208 nmax = std::max(nmax, nelm);
209 }
210 fXsec.resize(nmax);
211
212 if(isMaster) {
213 if(nullptr == theHandler) {
214 theHandler = new G4HadDataHandler(nTables);
215 }
216
217 fMaxEnergy = std::max(fMaxEnergy, param->GetMaxEnergy());
218 nLowE *= G4lrint(std::log10(fMiddleEnergy/fMinEnergy));
219 nHighE *= G4lrint(std::log10(fMaxEnergy/fMiddleEnergy));
220
221 G4PhysicsVector* vec = nullptr;
222 G4PhysicsLogVector aVector(fMinEnergy, fMiddleEnergy, nLowE, false);
223 G4PhysicsLogVector bVector(fMiddleEnergy, fMaxEnergy, nHighE, false);
224
225 for(std::size_t i=0; i<nTables; ++i) {
226 G4PhysicsTable* table = new G4PhysicsTable();
227 theHandler->UpdateTable(table, i);
228 table->resize(nmat);
229 for(std::size_t j=0; j<nmat; ++j) {
230 vec = (*table)[j];
231 if (nullptr == vec) {
232 if(i <= 2) {
233 vec = new G4PhysicsVector(aVector);
234 } else {
235 vec = new G4PhysicsVector(bVector);
236 }
238 }
239 }
240 }
241 }
242}
243
244//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
245
247{
248 if(1 < verboseLevel) {
249 G4cout << "### G4NeutronGeneralProcess::BuildPhysicsTable() for "
250 << GetProcessName()
251 << " and particle " << part.GetParticleName()
252 << G4endl;
253 }
254 fElastic->BuildPhysicsTable(part);
255 fInelastic->BuildPhysicsTable(part);
256 fCapture->BuildPhysicsTable(part);
257
258 if(isMaster) {
259 std::size_t nmat = G4Material::GetNumberOfMaterials();
261
262 auto tables = theHandler->GetTables();
263
264 G4double sigEl(0.), sigInel(0.), sigCap(0.), val(0.), sum(0.);
265
266 for(std::size_t i=0; i<nmat; ++i) {
267 const G4Material* mat = (*matTable)[i];
268
269 // energy interval 0
270 std::size_t nn = (*(tables[0]))[i]->GetVectorLength();
271 if(1 < verboseLevel) {
272 G4cout << "======= Zone 0 ======= N= " << nn
273 << " for " << mat->GetName() << G4endl;
274 }
275 for(std::size_t j=0; j<nn; ++j) {
276 G4double e = (*(tables[0]))[i]->Energy(j);
277 G4double loge = G4Log(e);
278 sigEl = fXSFactorEl*ComputeCrossSection(fElasticXS, mat, e, loge);
279 sigInel = fXSFactorInel*ComputeCrossSection(fInelasticXS, mat, e, loge);
280 sigCap = ComputeCrossSection(fCaptureXS, mat, e, loge);
281 sum = sigEl + sigInel + sigCap;
282 if(1 < verboseLevel) {
283 G4cout << j << ". E= " << e << " xs=" << sum << " sigEl=" << sigEl
284 << " sigInel=" << sigInel << " sigCap=" << sigCap << G4endl;
285 }
286 (*(tables[0]))[i]->PutValue(j, sum);
287 val = sigEl/sum;
288 (*(tables[1]))[i]->PutValue(j, val);
289 val = (sigEl + sigInel)/sum;
290 (*(tables[2]))[i]->PutValue(j, val);
291 }
292
293 // energy interval 1
294 nn = (*(tables[3]))[0]->GetVectorLength();
295 if(1 < verboseLevel) {
296 G4cout << "======= Zone 1 ======= N= " << nn << G4endl;
297 }
298 for(std::size_t j=0; j<nn; ++j) {
299 G4double e = (*(tables[3]))[i]->Energy(j);
300 G4double loge = G4Log(e);
301 sigEl = fXSFactorEl*ComputeCrossSection(fElasticXS, mat, e, loge);
302 sigInel = fXSFactorInel*ComputeCrossSection(fInelasticXS, mat, e, loge);
303 sum = sigEl + sigInel;
304 if(1 < verboseLevel) {
305 G4cout << j << ". E= " << e << " xs=" << sum << " sigEl=" << sigEl
306 << " sigInel=" << sigInel << " factInel=" << fXSFactorInel
307 << G4endl;
308 }
309 (*(tables[3]))[i]->PutValue(j, sum);
310 val = sigInel/sum;
311 (*(tables[4]))[i]->PutValue(j, val);
312 }
313 }
314 }
315 if(1 < verboseLevel) {
316 G4cout << "### G4VEmProcess::BuildPhysicsTable() done for "
317 << GetProcessName()
318 << " and particle " << part.GetParticleName()
319 << G4endl;
320 }
321}
322
323//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
324
326G4NeutronGeneralProcess::ComputeCrossSection(G4VCrossSectionDataSet* xs,
327 const G4Material* mat,
328 G4double e, G4double loge)
329{
330 const G4double* natom = mat->GetVecNbOfAtomsPerVolume();
331 G4int nelm = (G4int)mat->GetNumberOfElements();
332 G4double sig = 0.0;
333 for(G4int i=0; i<nelm; ++i) {
334 sig += natom[i]*xs->ComputeCrossSectionPerElement(e, loge, fNeutron,
335 mat->GetElement(i), mat);
336 }
337 return sig;
338}
339
340//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
341
343{
345 fCurrMat = nullptr;
346}
347
348//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
349
351 const G4Track& track,
352 G4double previousStepSize,
354{
356
357 // time limit
358 if(track.GetGlobalTime() >= fTimeLimit) {
359 fLambda = 0.0;
360 return 0.0;
361 }
362
363 // recompute total cross section if needed
364 CurrentCrossSection(track);
365
367
368 // beggining of tracking (or just after DoIt of this process)
371
372 } else {
373
375 previousStepSize/currentInteractionLength;
378 }
379
381 /*
382 G4cout << "PostStepGetPhysicalInteractionLength: e= " << energy
383 << " idxe= " << idxEnergy << " xs= " << fLambda
384 << " x= " << x << G4endl;
385 */
386 return x;
387}
388
389//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
390
392 const G4Step& step)
393{
394 fSelectedProc = this;
395 // time limit
396 if(0.0 == fLambda) {
399 return theTotalResult;
400 }
401 // In all cases clear number of interaction lengths
404 /*
405 G4cout << "PostStep: preStepLambda= " << fLambda << " idxE= " << idxEnergy
406 << " matIndex=" << matIndex << G4endl;
407 */
408 if (0 == idxEnergy) {
409 if(q <= GetProbability(1)) {
410 SelectedProcess(step, fElastic, fXSSElastic);
411 } else if(q <= GetProbability(2)) {
412 SelectedProcess(step, fInelastic, fXSSInelastic);
413 } else {
414 SelectedProcess(step, fCapture, fXSSCapture);
415 }
416 } else {
417 if(q <= GetProbability(4)) {
418 SelectedProcess(step, fInelastic, fXSSInelastic);
419 } else {
420 SelectedProcess(step, fElastic, fXSSElastic);
421 }
422 }
423 // total cross section is needed for selection of an element
424 if(fCurrMat->GetNumberOfElements() > 1) {
425 fCurrentXSS->ComputeCrossSection(track.GetDynamicParticle(), fCurrMat);
426 }
427 /*
428 G4cout << "## neutron E(MeV)=" << fCurrE << " inside " << fCurrMat->GetName()
429 << fSelectedProc->GetProcessName()
430 << " time(ns)=" << track.GetGlobalTime()/ns << G4endl;
431 */
432 // sample secondaries
433 return fSelectedProc->PostStepDoIt(track, step);
434}
435
436//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
437
438G4bool
440 const G4String& directory,
441 G4bool ascii)
442{
443 G4bool yes = true;
444 if(!isMaster) { return yes; }
445 for(std::size_t i=0; i<nTables; ++i) {
446 G4String nam = (0==i || 3==i)
447 ? "LambdaNeutronGeneral" + nameT[i] : "ProbNeutronGeneral" + nameT[i];
448 G4String fnam = GetPhysicsTableFileName(part, directory, nam, ascii);
449 auto table = theHandler->Table(i);
450 if(nullptr == table || !table->StorePhysicsTable(fnam, ascii)) {
451 yes = false;
452 }
453 }
454 return yes;
455}
456
457//....Ooooo0ooooo ........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
458
460 G4double,
462{
464 // recompute total cross section if needed
465 CurrentCrossSection(track);
467}
468
469//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
470
472{
473 fElastic->ProcessDescription(out);
474 fInelastic->ProcessDescription(out);
475 fCapture->ProcessDescription(out);
476}
477
478//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
479
481{
482 return (nullptr != fSelectedProc) ? fSelectedProc->GetProcessName()
484}
485
486//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
487
489{
490 return (nullptr != fSelectedProc) ? fSelectedProc->GetProcessSubType()
492}
493
494//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
495
497{
498 return fSelectedProc;
499}
500
501//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double condition(const G4ErrorSymMatrix &m)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ NotForced
@ fNeutronGeneral
G4double G4Log(G4double x)
Definition: G4Log.hh:227
std::vector< G4Material * > G4MaterialTable
@ 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
const std::vector< G4VCrossSectionDataSet * > & GetDataSetList() const
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
G4PhysicsTable * Table(std::size_t idx) const
void UpdateTable(G4PhysicsTable *, std::size_t idx)
const std::vector< G4PhysicsTable * > & GetTables() const
static G4HadronicParameters * Instance()
G4double XSFactorNucleonElastic() const
G4double GetMaxEnergy() const
G4double XSFactorNucleonInelastic() const
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
void ProcessDescription(std::ostream &outFile) const override
G4ParticleChange * theTotalResult
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
void PreparePhysicsTable(const G4ParticleDefinition &) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4CrossSectionDataStore * GetCrossSectionDataStore()
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:684
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:197
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:201
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677
const G4String & GetName() const
Definition: G4Material.hh:172
void SetCaptureProcess(G4HadronicProcess *)
G4NeutronGeneralProcess(const G4String &pname="NeutronGeneralProc")
void StartTracking(G4Track *) override
void ProcessDescription(std::ostream &outFile) const override
void PreparePhysicsTable(const G4ParticleDefinition &) override
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
void SelectedProcess(const G4Step &step, G4HadronicProcess *ptr, G4CrossSectionDataStore *)
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
void SetElasticProcess(G4HadronicProcess *)
G4double GetProbability(size_t idxt)
G4bool StorePhysicsTable(const G4ParticleDefinition *part, const G4String &directory, G4bool ascii) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
const G4VProcess * GetCreatorProcess() const override
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
const G4String & GetSubProcessName() const
void SetInelasticProcess(G4HadronicProcess *)
G4bool IsApplicable(const G4ParticleDefinition &) override
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
void Initialize(const G4Track &) override
const G4String & GetParticleName() const
static void SetPhysicsVector(G4PhysicsTable *physTable, std::size_t idx, G4PhysicsVector *vec)
void resize(std::size_t, G4PhysicsVector *vec=nullptr)
Definition: G4Step.hh:62
G4double GetGlobalTime() const
const G4DynamicParticle * GetDynamicParticle() const
virtual G4double ComputeCrossSectionPerElement(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *mat=nullptr)
void ProposeTrackStatus(G4TrackStatus status)
G4double currentInteractionLength
Definition: G4VProcess.hh:339
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:342
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:416
G4int verboseLevel
Definition: G4VProcess.hh:360
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:335
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
G4int GetProcessSubType() const
Definition: G4VProcess.hh:404
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:188
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
Definition: DoubConv.h:17
G4bool IsWorkerThread()
Definition: G4Threading.cc:123
int G4lrint(double ad)
Definition: templates.hh:134