source: trunk/source/processes/hadronic/management/src/G4HadronicProcessStore.cc@ 1330

Last change on this file since 1330 was 1315, checked in by garnier, 15 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 21.3 KB
Line 
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: G4HadronicProcessStore.cc,v 1.14 2010/04/21 17:55:25 dennis Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-cand-01 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4HadronicProcessStore
35//
36// Author: Vladimir Ivanchenko
37//
38// Creation date: 09.05.2008
39//
40// Modifications:
41// 23.01.2009 V.Ivanchenko add destruction of processes
42//
43// Class Description:
44// Singleton to store hadronic processes, to provide access to processes
45// and to printout information about processes
46//
47// -------------------------------------------------------------------
48//
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
51
52#include "G4HadronicProcessStore.hh"
53#include "G4Element.hh"
54#include "G4ProcessManager.hh"
55#include "G4Electron.hh"
56#include "G4Proton.hh"
57#include "G4HadronicInteractionRegistry.hh"
58#include "G4CrossSectionDataSetRegistry.hh"
59#include "G4HadronicEPTestMessenger.hh"
60
61G4HadronicProcessStore* G4HadronicProcessStore::theInstance = 0;
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
64
65G4HadronicProcessStore* G4HadronicProcessStore::Instance()
66{
67 if(0 == theInstance) {
68 static G4HadronicProcessStore manager;
69 theInstance = &manager;
70 }
71 return theInstance;
72}
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
75
76G4HadronicProcessStore::~G4HadronicProcessStore()
77{
78 Clean();
79 G4HadronicInteractionRegistry::Instance()->Clean();
80 G4CrossSectionDataSetRegistry::Instance()->Clean();
81 delete theEPTestMessenger;
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
85
86void G4HadronicProcessStore::Clean()
87{
88 G4int i;
89 //G4cout << "G4HadronicProcessStore::Clean() Nproc= " << n_proc
90 // << " Nextra= " << n_extra << G4endl;
91 if(n_proc > 0) {
92 for (i=0; i<n_proc; i++) {
93 if( process[i] ) {
94 //G4cout << "G4HadronicProcessStore::Clean() delete hadronic " << i << G4endl;
95 //G4cout << process[i]->GetProcessName() << G4endl;
96 G4HadronicProcess* p = process[i];
97 process[i] = 0;
98 delete p;
99 }
100 }
101 }
102 if(n_extra > 0) {
103 for(i=0; i<n_extra; i++) {
104 if(extraProcess[i]) {
105 //G4cout << "G4HadronicProcessStore::Clean() delete extra "
106 // << i << G4endl;
107 //G4cout << extraProcess[i]->GetProcessName() << G4endl;
108 G4VProcess* p = extraProcess[i];
109 extraProcess[i] = 0;
110 delete p;
111 }
112 }
113 }
114 //G4cout << "G4HadronicProcessStore::Clean() done" << G4endl;
115 n_extra = 0;
116 n_proc = 0;
117}
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
120
121G4HadronicProcessStore::G4HadronicProcessStore()
122{
123 n_proc = 0;
124 n_part = 0;
125 n_model= 0;
126 n_extra= 0;
127 currentProcess = 0;
128 currentParticle = 0;
129 verbose = 1;
130 buildTableStart = true;
131 theEPTestMessenger = new G4HadronicEPTestMessenger(this);
132}
133
134//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
135
136G4double G4HadronicProcessStore::GetElasticCrossSectionPerVolume(
137 const G4ParticleDefinition *aParticle,
138 G4double kineticEnergy,
139 const G4Material *material)
140{
141 G4double cross = 0.0;
142 const G4ElementVector* theElementVector = material->GetElementVector();
143 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
144 size_t nelm = material->GetNumberOfElements();
145 for (size_t i=0; i<nelm; i++) {
146 const G4Element* elm = (*theElementVector)[i];
147 cross += theAtomNumDensityVector[i]*
148 GetElasticCrossSectionPerAtom(aParticle,kineticEnergy,elm);
149 }
150 return cross;
151}
152
153//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
154
155G4double G4HadronicProcessStore::GetElasticCrossSectionPerAtom(
156 const G4ParticleDefinition *aParticle,
157 G4double kineticEnergy,
158 const G4Element *anElement)
159{
160 G4HadronicProcess* hp = FindProcess(aParticle, fHadronElastic);
161 localDP.SetKineticEnergy(kineticEnergy);
162 G4double cross = 0.0;
163 if(hp) cross = hp->GetMicroscopicCrossSection(&localDP,
164 anElement,
165 STP_Temperature);
166 return cross;
167}
168
169//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
170
171G4double G4HadronicProcessStore::GetElasticCrossSectionPerIsotope(
172 const G4ParticleDefinition*,
173 G4double,
174 G4int, G4int)
175{
176 return 0.0;
177}
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
180
181G4double G4HadronicProcessStore::GetInelasticCrossSectionPerVolume(
182 const G4ParticleDefinition *aParticle,
183 G4double kineticEnergy,
184 const G4Material *material)
185{
186 G4double cross = 0.0;
187 const G4ElementVector* theElementVector = material->GetElementVector();
188 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
189 size_t nelm = material->GetNumberOfElements();
190 for (size_t i=0; i<nelm; i++) {
191 const G4Element* elm = (*theElementVector)[i];
192 cross += theAtomNumDensityVector[i]*
193 GetInelasticCrossSectionPerAtom(aParticle,kineticEnergy,elm);
194 }
195 return cross;
196}
197
198//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
199
200G4double G4HadronicProcessStore::GetInelasticCrossSectionPerAtom(
201 const G4ParticleDefinition *aParticle,
202 G4double kineticEnergy,
203 const G4Element *anElement)
204{
205 G4HadronicProcess* hp = FindProcess(aParticle, fHadronInelastic);
206 localDP.SetDefinition(const_cast<G4ParticleDefinition*>(aParticle));
207 localDP.SetKineticEnergy(kineticEnergy);
208 G4double cross = 0.0;
209 if(hp) cross = hp->GetMicroscopicCrossSection(&localDP,
210 anElement,
211 STP_Temperature);
212 return cross;
213}
214
215//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
216
217G4double G4HadronicProcessStore::GetInelasticCrossSectionPerIsotope(
218 const G4ParticleDefinition *,
219 G4double,
220 G4int, G4int)
221{
222 return 0.0;
223}
224
225//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
226
227G4double G4HadronicProcessStore::GetCaptureCrossSectionPerVolume(
228 const G4ParticleDefinition *aParticle,
229 G4double kineticEnergy,
230 const G4Material *material)
231{
232 G4double cross = 0.0;
233 const G4ElementVector* theElementVector = material->GetElementVector();
234 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
235 size_t nelm = material->GetNumberOfElements();
236 for (size_t i=0; i<nelm; i++) {
237 const G4Element* elm = (*theElementVector)[i];
238 cross += theAtomNumDensityVector[i]*
239 GetCaptureCrossSectionPerAtom(aParticle,kineticEnergy,elm);
240 }
241 return cross;
242}
243
244//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
245
246G4double G4HadronicProcessStore::GetCaptureCrossSectionPerAtom(
247 const G4ParticleDefinition *aParticle,
248 G4double kineticEnergy,
249 const G4Element *anElement)
250{
251 G4HadronicProcess* hp = FindProcess(aParticle, fCapture);
252 localDP.SetDefinition(const_cast<G4ParticleDefinition*>(aParticle));
253 localDP.SetKineticEnergy(kineticEnergy);
254 G4double cross = 0.0;
255 if(hp) cross = hp->GetMicroscopicCrossSection(&localDP,
256 anElement,
257 STP_Temperature);
258 return cross;
259}
260
261//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
262
263G4double G4HadronicProcessStore::GetCaptureCrossSectionPerIsotope(
264 const G4ParticleDefinition *,
265 G4double,
266 G4int, G4int)
267{
268 return 0.0;
269}
270
271//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
272
273G4double G4HadronicProcessStore::GetFissionCrossSectionPerVolume(
274 const G4ParticleDefinition *aParticle,
275 G4double kineticEnergy,
276 const G4Material *material)
277{
278 G4double cross = 0.0;
279 const G4ElementVector* theElementVector = material->GetElementVector();
280 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
281 size_t nelm = material->GetNumberOfElements();
282 for (size_t i=0; i<nelm; i++) {
283 const G4Element* elm = (*theElementVector)[i];
284 cross += theAtomNumDensityVector[i]*
285 GetFissionCrossSectionPerAtom(aParticle,kineticEnergy,elm);
286 }
287 return cross;
288}
289
290//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
291
292G4double G4HadronicProcessStore::GetFissionCrossSectionPerAtom(
293 const G4ParticleDefinition *aParticle,
294 G4double kineticEnergy,
295 const G4Element *anElement)
296{
297 G4HadronicProcess* hp = FindProcess(aParticle, fFission);
298 localDP.SetDefinition(const_cast<G4ParticleDefinition*>(aParticle));
299 localDP.SetKineticEnergy(kineticEnergy);
300 G4double cross = 0.0;
301 if(hp) cross = hp->GetMicroscopicCrossSection(&localDP,
302 anElement,
303 STP_Temperature);
304 return cross;
305}
306
307//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
308
309G4double G4HadronicProcessStore::GetFissionCrossSectionPerIsotope(
310 const G4ParticleDefinition *,
311 G4double,
312 G4int, G4int)
313{
314 return 0.0;
315}
316
317//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
318
319G4double G4HadronicProcessStore::GetChargeExchangeCrossSectionPerVolume(
320 const G4ParticleDefinition *aParticle,
321 G4double kineticEnergy,
322 const G4Material *material)
323{
324 G4double cross = 0.0;
325 const G4ElementVector* theElementVector = material->GetElementVector();
326 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
327 size_t nelm = material->GetNumberOfElements();
328 for (size_t i=0; i<nelm; i++) {
329 const G4Element* elm = (*theElementVector)[i];
330 cross += theAtomNumDensityVector[i]*
331 GetChargeExchangeCrossSectionPerAtom(aParticle,kineticEnergy,elm);
332 }
333 return cross;
334}
335
336//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
337
338G4double G4HadronicProcessStore::GetChargeExchangeCrossSectionPerAtom(
339 const G4ParticleDefinition *aParticle,
340 G4double kineticEnergy,
341 const G4Element *anElement)
342{
343 G4HadronicProcess* hp = FindProcess(aParticle, fChargeExchange);
344 localDP.SetDefinition(const_cast<G4ParticleDefinition*>(aParticle));
345 localDP.SetKineticEnergy(kineticEnergy);
346 G4double cross = 0.0;
347 if(hp) cross = hp->GetMicroscopicCrossSection(&localDP,
348 anElement,
349 STP_Temperature);
350 return cross;
351}
352
353//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
354
355G4double G4HadronicProcessStore::GetChargeExchangeCrossSectionPerIsotope(
356 const G4ParticleDefinition *,
357 G4double,
358 G4int, G4int)
359{
360 return 0.0;
361}
362
363//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
364
365void G4HadronicProcessStore::Register(G4HadronicProcess* proc)
366{
367 if(0 < n_proc) {
368 for(G4int i=0; i<n_proc; i++) {
369 if(process[i] == proc) { return; }
370 }
371 }
372 // G4cout << "G4HadronicProcessStore::Register hadronic " << n_proc
373 // << " " << proc->GetProcessName() << G4endl;
374 n_proc++;
375 process.push_back(proc);
376}
377
378//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
379
380void G4HadronicProcessStore::RegisterParticle(G4HadronicProcess* proc,
381 const G4ParticleDefinition* part)
382{
383 G4int i=0;
384 for(; i<n_proc; i++) {if(process[i] == proc) break;}
385 G4int j=0;
386 for(; j<n_part; j++) {if(particle[j] == part) break;}
387
388 if(j == n_part) {
389 n_part++;
390 particle.push_back(part);
391 wasPrinted.push_back(0);
392 }
393
394 // the pair should be added?
395 if(i < n_proc) {
396 std::multimap<PD,HP,std::less<PD> >::iterator it;
397 for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
398 if(it->first == part) {
399 HP process = (it->second);
400 if(proc == process) return;
401 }
402 }
403 }
404
405 p_map.insert(std::multimap<PD,HP>::value_type(part,proc));
406}
407
408//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
409
410void G4HadronicProcessStore::RegisterInteraction(G4HadronicProcess* proc,
411 G4HadronicInteraction* mod)
412{
413 G4int i=0;
414 for(; i<n_proc; i++) {if(process[i] == proc) break;}
415 G4int k=0;
416 for(; k<n_model; k++) {if(model[k] == mod) break;}
417
418 m_map.insert(std::multimap<HP,HI>::value_type(proc,mod));
419
420 if(k == n_model) {
421 n_model++;
422 model.push_back(mod);
423 modelName.push_back(mod->GetModelName());
424 }
425}
426
427//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
428
429void G4HadronicProcessStore::DeRegister(G4HadronicProcess* proc)
430{
431 if(0 == n_proc) return;
432 for(G4int i=0; i<n_proc; i++) {
433 if(process[i] == proc) {
434 process[i] = 0;
435 return;
436 }
437 }
438}
439
440//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
441
442void G4HadronicProcessStore::RegisterExtraProcess(G4VProcess* proc)
443{
444 if(0 < n_extra) {
445 for(G4int i=0; i<n_extra; i++) {
446 if(extraProcess[i] == proc) return;
447 }
448 }
449 //G4cout << "Extra Process: " << n_extra << " " << proc->GetProcessName()
450 // << " " << proc << G4endl;
451
452 n_extra++;
453 extraProcess.push_back(proc);
454}
455
456//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
457
458void G4HadronicProcessStore::RegisterParticleForExtraProcess(
459 G4VProcess* proc,
460 const G4ParticleDefinition* part)
461{
462 G4int i=0;
463 for(; i<n_extra; i++) {if(extraProcess[i] == proc) break;}
464 G4int j=0;
465 for(; j<n_part; j++) {if(particle[j] == part) break;}
466
467 if(j == n_part) {
468 n_part++;
469 particle.push_back(part);
470 wasPrinted.push_back(0);
471 }
472
473 // the pair should be added?
474 if(i < n_extra) {
475 std::multimap<PD,G4VProcess*,std::less<PD> >::iterator it;
476 for(it=ep_map.lower_bound(part); it!=ep_map.upper_bound(part); ++it) {
477 if(it->first == part) {
478 G4VProcess* process = (it->second);
479 if(proc == process) return;
480 }
481 }
482 }
483
484 ep_map.insert(std::multimap<PD,G4VProcess*>::value_type(part,proc));
485}
486
487//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
488
489void G4HadronicProcessStore::DeRegisterExtraProcess(G4VProcess* proc)
490{
491 //G4cout << "Deregister Extra Process: " << proc << " " << proc->GetProcessName() << G4endl;
492 if(0 == n_extra) return;
493 for(G4int i=0; i<n_extra; i++) {
494 if(extraProcess[i] == proc) {
495 extraProcess[i] = 0;
496 //G4cout << "Extra Process: " << i << " is deregisted " << G4endl;
497 return;
498 }
499 }
500}
501
502//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
503
504void G4HadronicProcessStore::PrintInfo(const G4ParticleDefinition* part)
505{
506 if(buildTableStart && part == particle[n_part - 1]) {
507 buildTableStart = false;
508 Dump(verbose);
509 }
510}
511
512//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
513
514void G4HadronicProcessStore::Dump(G4int level)
515{
516 if(level > 0) {
517 G4cout << "=============================================================="
518 << "=============================="
519 << G4endl;
520 G4cout << " HADRONIC PROCESSES SUMMARY (verbose level " << level
521 << ")" << G4endl;
522 }
523 for(G4int i=0; i<n_part; i++) {
524 PD part = particle[i];
525 G4String pname = part->GetParticleName();
526 G4bool yes = false;
527 if(level >= 2) yes = true;
528 else if(level == 1 && (pname == "proton" ||
529 pname == "neutron" ||
530 pname == "pi+" ||
531 pname == "pi-" ||
532 pname == "gamma" ||
533 pname == "e-" ||
534 pname == "mu-" ||
535 pname == "kaon+" ||
536 pname == "kaon-" ||
537 pname == "lambda" ||
538 pname == "anti_neutron" ||
539 pname == "anti_proton")) yes = true;
540 if(yes) {
541 // main processes
542 std::multimap<PD,HP,std::less<PD> >::iterator it;
543 for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
544 if(it->first == part) {
545 HP proc = (it->second);
546 G4int j=0;
547 for(; j<n_proc; j++) {
548 if(process[j] == proc) {
549 Print(j, i);
550 }
551 }
552 }
553 }
554 // extra processes
555 std::multimap<PD,G4VProcess*,std::less<PD> >::iterator itp;
556 for(itp=ep_map.lower_bound(part); itp!=ep_map.upper_bound(part); ++itp) {
557 if(itp->first == part) {
558 G4VProcess* proc = (itp->second);
559 if(wasPrinted[i] == 0) {
560 wasPrinted[i] = 1;
561 G4cout<<G4endl;
562 G4cout << " Hadronic Processes for <"
563 <<part->GetParticleName() << ">" << G4endl;
564 }
565 G4cout << " " << proc->GetProcessName() << G4endl;
566 }
567 }
568 }
569 }
570 if(level > 0) {
571 G4cout << "=============================================================="
572 << "=============================="
573 << G4endl;
574 }
575}
576
577//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
578
579void G4HadronicProcessStore::Print(G4int idxProc, G4int idxPart)
580{
581 G4HadronicProcess* proc = process[idxProc];
582 const G4ParticleDefinition* part = particle[idxPart];
583 if(wasPrinted[idxPart] == 0) {
584 wasPrinted[idxPart] = 1;
585 G4cout<<G4endl;
586 G4cout << " Hadronic Processes for <"
587 <<part->GetParticleName() << ">" << G4endl;
588 }
589 HI hi = 0;
590 G4bool first;
591 std::multimap<HP,HI,std::less<HP> >::iterator ih;
592 G4cout << std::setw(20) << proc->GetProcessName()
593 << " Models: ";
594 first = true;
595 for(ih=m_map.lower_bound(proc); ih!=m_map.upper_bound(proc); ++ih) {
596 if(ih->first == proc) {
597 hi = ih->second;
598 G4int i=0;
599 for(; i<n_model; i++) {
600 if(model[i] == hi) break;
601 }
602 if(!first) G4cout << " ";
603 first = false;
604 G4cout << std::setw(25) << modelName[i]
605 << ": Emin(GeV)= "
606 << std::setw(5) << hi->GetMinEnergy()/GeV
607 << " Emax(GeV)= "
608 << hi->GetMaxEnergy()/GeV
609 << G4endl;
610 }
611 }
612}
613
614//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
615
616void G4HadronicProcessStore::SetVerbose(G4int val)
617{
618 verbose = val;
619 G4int i;
620 for(i=0; i<n_proc; i++) {
621 if(process[i]) process[i]->SetVerboseLevel(val);
622 }
623 for(i=0; i<n_model; i++) {
624 if(model[i]) model[i]->SetVerboseLevel(val);
625 }
626}
627
628//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
629
630G4int G4HadronicProcessStore::GetVerbose()
631{
632 return verbose;
633}
634
635//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
636
637G4HadronicProcess* G4HadronicProcessStore::FindProcess(
638 const G4ParticleDefinition* part, G4HadronicProcessType subType)
639{
640 bool isNew = false;
641 G4HadronicProcess* hp = 0;
642
643 if(part != currentParticle) {
644 isNew = true;
645 currentParticle = part;
646 localDP.SetDefinition(const_cast<G4ParticleDefinition*>(part));
647 } else if(!currentProcess) {
648 isNew = true;
649 } else if(subType == currentProcess->GetProcessSubType()) {
650 hp = currentProcess;
651 } else {
652 isNew = true;
653 }
654
655 if(isNew) {
656 std::multimap<PD,HP,std::less<PD> >::iterator it;
657 for(it=p_map.lower_bound(part); it!=p_map.upper_bound(part); ++it) {
658 if(it->first == part && subType == (it->second)->GetProcessSubType()) {
659 hp = it->second;
660 break;
661 }
662 }
663 currentProcess = hp;
664 }
665
666 return hp;
667}
668
669
670void G4HadronicProcessStore::SetEpReportLevel(G4int level)
671{
672 G4cout << " Setting energy/momentum report level to " << level
673 << " for " << process.size() << " hadronic processes " << G4endl;
674 for (G4int i = 0; i < G4int(process.size()); i++) {
675 process[i]->SetEpReportLevel(level);
676 }
677}
678
679
680void G4HadronicProcessStore::SetProcessAbsLevel(G4double abslevel)
681{
682 G4cout << " Setting absolute energy/momentum test level to " << abslevel << G4endl;
683 G4double rellevel = 0.0;
684 G4HadronicProcess* theProcess = 0;
685 for (G4int i = 0; i < G4int(process.size()); i++) {
686 theProcess = process[i];
687 rellevel = theProcess->GetEnergyMomentumCheckLevels().first;
688 theProcess->SetEnergyMomentumCheckLevels(rellevel, abslevel);
689 }
690}
691
692
693void G4HadronicProcessStore::SetProcessRelLevel(G4double rellevel)
694{
695 G4cout << " Setting relative energy/momentum test level to " << rellevel << G4endl;
696 G4double abslevel = 0.0;
697 G4HadronicProcess* theProcess = 0;
698 for (G4int i = 0; i < G4int(process.size()); i++) {
699 theProcess = process[i];
700 abslevel = theProcess->GetEnergyMomentumCheckLevels().second;
701 theProcess->SetEnergyMomentumCheckLevels(rellevel, abslevel);
702 }
703}
Note: See TracBrowser for help on using the repository browser.