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

Last change on this file since 1228 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

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