Loading...
Searching...
No Matches
B03PhysicsList.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/// \file biasing/B03/src/B03PhysicsList.cc
27/// \brief Implementation of the B03PhysicsList class
28//
29//
30//
31
32#include "globals.hh"
33#include <iomanip>
34
35#include "B03PhysicsList.hh"
36
37#include "G4ParticleDefinition.hh"
38#include "G4ParticleWithCuts.hh"
39#include "G4ProcessManager.hh"
40#include "G4ProcessVector.hh"
41#include "G4ParticleTypes.hh"
42#include "G4ParticleTable.hh"
43#include "G4BosonConstructor.hh"
44#include "G4LeptonConstructor.hh"
45#include "G4MesonConstructor.hh"
46#include "G4BaryonConstructor.hh"
47#include "G4IonConstructor.hh"
48#include "G4ShortLivedConstructor.hh"
49#include "G4Material.hh"
50#include "G4MaterialTable.hh"
51#include "G4SystemOfUnits.hh"
52#include "G4HadronicParameters.hh"
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55
58fBiasWorldName(parallelname)
59{
60 fParaWorldName.clear();
61 SetVerboseLevel(1);
62}
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65
70
71//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72
74{
75 // In this method, static member functions should be called
76 // for all particles which you want to use.
77 // This ensures that objects of these particle types will be
78 // created in the program.
79
86}
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
89
91{
92 // Construct all bosons
93 G4BosonConstructor pConstructor;
94 pConstructor.ConstructParticle();
95}
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
98
100{
101 // Construct all leptons
102 G4LeptonConstructor pConstructor;
103 pConstructor.ConstructParticle();
104}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107
109{
110 // Construct all mesons
111 G4MesonConstructor pConstructor;
112 pConstructor.ConstructParticle();
113}
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116
118{
119 // Construct all barions
120 G4BaryonConstructor pConstructor;
121 pConstructor.ConstructParticle();
122}
123
124//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125
127{
128 // Construct light ions
129 G4IonConstructor pConstructor;
130 pConstructor.ConstructParticle();
131}
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134
136{
137 // Construct resonaces and quarks
138 G4ShortLivedConstructor pConstructor;
139 pConstructor.ConstructParticle();
140}
141
142//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
143
145{
146 AddTransportation();
149 ConstructEM();
151 ConstructHad();
153}
154
155//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
156
157#include "G4ComptonScattering.hh"
158#include "G4GammaConversion.hh"
159#include "G4PhotoElectricEffect.hh"
160
161#include "G4eMultipleScattering.hh"
162#include "G4MuMultipleScattering.hh"
163#include "G4hMultipleScattering.hh"
164
165#include "G4eIonisation.hh"
166#include "G4eBremsstrahlung.hh"
167#include "G4eplusAnnihilation.hh"
168
169#include "G4MuIonisation.hh"
170#include "G4MuBremsstrahlung.hh"
171#include "G4MuPairProduction.hh"
172
173#include "G4hIonisation.hh"
174
176{
177 auto particleIterator=GetParticleIterator();
178 particleIterator->reset();
179 while( (*particleIterator)() ){
180 G4ParticleDefinition* particle = particleIterator->value();
181 G4ProcessManager* pmanager = particle->GetProcessManager();
182 G4String particleName = particle->GetParticleName();
183
184 if (particleName == "gamma") {
185 // gamma
186 // Construct processes for gamma
187 pmanager->AddDiscreteProcess(new G4GammaConversion());
188 pmanager->AddDiscreteProcess(new G4ComptonScattering());
189 pmanager->AddDiscreteProcess(new G4PhotoElectricEffect());
190
191 } else if (particleName == "e-") {
192 //electron
193 // Construct processes for electron
194 pmanager->AddProcess(new G4eMultipleScattering(),-1,1,1);
195 pmanager->AddProcess(new G4eIonisation(),-1,2,2);
196 pmanager->AddProcess(new G4eBremsstrahlung(),-1,-1,3);
197
198 } else if (particleName == "e+") {
199 //positron
200 // Construct processes for positron
201 pmanager->AddProcess(new G4eMultipleScattering(),-1,1,1);
202
203 pmanager->AddProcess(new G4eIonisation(),-1,2,2);
204 pmanager->AddProcess(new G4eBremsstrahlung(),-1,-1,3);
205 pmanager->AddProcess(new G4eplusAnnihilation(),0,-1,4);
206
207 } else if( particleName == "mu+" ||
208 particleName == "mu-" ) {
209 //muon
210 // Construct processes for muon+
211 pmanager->AddProcess(new G4MuMultipleScattering(),-1,1,1);
212 pmanager->AddProcess(new G4MuIonisation(),-1,2,2);
213 pmanager->AddProcess(new G4MuBremsstrahlung(),-1,-1,3);
214 pmanager->AddProcess(new G4MuPairProduction(),-1,-1,4);
215
216 } else if( particleName == "GenericIon" ) {
217 pmanager->AddProcess(new G4hMultipleScattering(),-1,1,1);
218 pmanager->AddProcess(new G4hIonisation(),-1,2,2);
219 } else {
220 if ((particle->GetPDGCharge() != 0.0) &&
221 (particle->GetParticleName() != "chargedgeantino")&&
222 (!particle->IsShortLived()) ) {
223 // all others charged particles except geantino
224 pmanager->AddProcess(new G4hMultipleScattering(),-1,1,1);
225 pmanager->AddProcess(new G4hIonisation(),-1,2,2);
226 }
227 }
228 }
229}
230
231//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
232
233// Hadron Processes
234
235#include "G4HadronElasticProcess.hh"
236#include "G4NeutronFissionProcess.hh"
237#include "G4NeutronCaptureProcess.hh"
238#include "G4HadronInelasticProcess.hh"
239
240// Low-energy Models
241
242#include "G4HadronElastic.hh"
243#include "G4NeutronRadCapture.hh"
244#include "G4LFission.hh"
245
246// -- generator models
247#include "G4TheoFSGenerator.hh"
248#include "G4ExcitationHandler.hh"
249#include "G4CompetitiveFission.hh"
250#include "G4GeneratorPrecompoundInterface.hh"
251#include "G4Fancy3DNucleus.hh"
252#include "G4CascadeInterface.hh"
253#include "G4StringModel.hh"
254#include "G4PreCompoundModel.hh"
255#include "G4FTFModel.hh"
256#include "G4QGSMFragmentation.hh"
257#include "G4LundStringFragmentation.hh"
258#include "G4ExcitedStringDecay.hh"
259#include "G4QMDReaction.hh"
260#include "G4BinaryLightIonReaction.hh"
261
262// Cross sections
263#include "G4ComponentGGNuclNuclXsc.hh"
264#include "G4ComponentGGHadronNucleusXsc.hh"
265#include "G4CrossSectionElastic.hh"
266#include "G4CrossSectionInelastic.hh"
267#include "G4NeutronInelasticXS.hh"
268
269//
270// ConstructHad()
271//
272// Makes discrete physics processes for the hadrons
273// The processes are: Elastic scattering, Inelastic scattering,
274// Fission (for neutron only), and Capture (neutron).
275//
276
278{
279 // this will be the model class for high energies
280 G4TheoFSGenerator* theTheoModel = new G4TheoFSGenerator;
281 G4TheoFSGenerator* antiBHighEnergyModel = new G4TheoFSGenerator;
282
283 // Evaporation logic
285 theHandler->SetMinEForMultiFrag(3*MeV);
286
287 // Pre equilibrium stage
288 G4PreCompoundModel* thePreEquilib = new G4PreCompoundModel(theHandler);
289
290 // a no-cascade generator-precompound interaface
293 theCascade->SetDeExcitation(thePreEquilib);
294
295 // Bertini cascade
297 bertini->SetMaxEnergy(22*MeV);
298
299 // here come the high energy parts
300 G4VPartonStringModel* theStringModel;
301 theStringModel = new G4FTFModel;
302 theTheoModel->SetTransport(theCascade);
303 theTheoModel->SetHighEnergyGenerator(theStringModel);
304 theTheoModel->SetMinEnergy(19*GeV);
305 theTheoModel->SetMaxEnergy( G4HadronicParameters::Instance()->GetMaxEnergy() );
306
307 G4VLongitudinalStringDecay* theFragmentation = new G4QGSMFragmentation;
308 G4ExcitedStringDecay* theStringDecay =
309 new G4ExcitedStringDecay(theFragmentation);
310 theStringModel->SetFragmentationModel(theStringDecay);
311
312 // high energy model for anti-baryons
313 antiBHighEnergyModel = new G4TheoFSGenerator("ANTI-FTFP");
314 G4FTFModel* antiBStringModel = new G4FTFModel;
315 G4ExcitedStringDecay* stringDecay =
317 antiBStringModel->SetFragmentationModel(stringDecay);
318
319 G4GeneratorPrecompoundInterface* antiBCascade =
321 G4PreCompoundModel* preEquilib =
323 antiBCascade->SetDeExcitation(preEquilib);
324
325 antiBHighEnergyModel->SetTransport(antiBCascade);
326 antiBHighEnergyModel->SetHighEnergyGenerator(antiBStringModel);
327 antiBHighEnergyModel->SetMinEnergy(0.0);
328 antiBHighEnergyModel->SetMaxEnergy(20*TeV);
329
330 // Light ion models
332 binaryCascade->SetMinEnergy(0.0);
333 binaryCascade->SetMaxEnergy(110*MeV);
334
335 G4QMDReaction* qmd = new G4QMDReaction;
336 qmd->SetMinEnergy(100*MeV);
337 qmd->SetMaxEnergy(10*GeV);
338
340
342 G4VCrossSectionDataSet * theGGHNEl = new G4CrossSectionElastic(ggHNXsec);
343 G4VCrossSectionDataSet * theGGHNInel = new G4CrossSectionInelastic(ggHNXsec);
344
345 // Elastic process
346 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
347 theElasticProcess->AddDataSet(theGGHNEl);
348 G4HadronElastic* theElasticModel = new G4HadronElastic;
349 theElasticProcess->RegisterMe(theElasticModel);
350
351 auto particleIterator=GetParticleIterator();
352 particleIterator->reset();
353 while ((*particleIterator)()) {
354 G4ParticleDefinition* particle = particleIterator->value();
355 G4ProcessManager* pmanager = particle->GetProcessManager();
356 G4String particleName = particle->GetParticleName();
357
358 if (particleName == "pi+") {
359 pmanager->AddDiscreteProcess(theElasticProcess);
360 G4HadronInelasticProcess* theInelasticProcess =
361 new G4HadronInelasticProcess( "inelastic", G4PionPlus::Definition() );
362 theInelasticProcess->AddDataSet(theGGHNInel);
363 theInelasticProcess->RegisterMe(bertini);
364 theInelasticProcess->RegisterMe(theTheoModel);
365 pmanager->AddDiscreteProcess(theInelasticProcess);
366 } else if (particleName == "pi-") {
367 pmanager->AddDiscreteProcess(theElasticProcess);
368 G4HadronInelasticProcess* theInelasticProcess =
369 new G4HadronInelasticProcess( "inelastic", G4PionMinus::Definition() );
370 theInelasticProcess->AddDataSet(theGGHNInel);
371 theInelasticProcess->RegisterMe(bertini);
372 theInelasticProcess->RegisterMe(theTheoModel);
373 pmanager->AddDiscreteProcess(theInelasticProcess);
374 } else if (particleName == "kaon+") {
375 pmanager->AddDiscreteProcess(theElasticProcess);
376 G4HadronInelasticProcess* theInelasticProcess =
377 new G4HadronInelasticProcess( "inelastic", G4KaonPlus::Definition() );
378 theInelasticProcess->AddDataSet(theGGHNInel);
379 theInelasticProcess->RegisterMe(bertini);
380 theInelasticProcess->RegisterMe(theTheoModel);
381 pmanager->AddDiscreteProcess(theInelasticProcess);
382 }
383 else if (particleName == "kaon0S") {
384 pmanager->AddDiscreteProcess(theElasticProcess);
385 G4HadronInelasticProcess* theInelasticProcess =
386 new G4HadronInelasticProcess( "inelastic", G4KaonZeroShort::Definition() );
387 theInelasticProcess->AddDataSet(theGGHNInel);
388 theInelasticProcess->RegisterMe(bertini);
389 theInelasticProcess->RegisterMe(theTheoModel);
390 pmanager->AddDiscreteProcess(theInelasticProcess);
391 }
392 else if (particleName == "kaon0L") {
393 pmanager->AddDiscreteProcess(theElasticProcess);
394 G4HadronInelasticProcess* theInelasticProcess =
395 new G4HadronInelasticProcess( "inelastic", G4KaonZeroLong::Definition() );
396 theInelasticProcess->AddDataSet(theGGHNInel);
397 theInelasticProcess->RegisterMe(bertini);
398 theInelasticProcess->RegisterMe(theTheoModel);
399 pmanager->AddDiscreteProcess(theInelasticProcess);
400 }
401 else if (particleName == "kaon-") {
402 pmanager->AddDiscreteProcess(theElasticProcess);
403 G4HadronInelasticProcess* theInelasticProcess =
404 new G4HadronInelasticProcess( "inelastic", G4KaonMinus::Definition() );
405 theInelasticProcess->AddDataSet(theGGHNInel);
406 theInelasticProcess->RegisterMe(bertini);
407 theInelasticProcess->RegisterMe(theTheoModel);
408 pmanager->AddDiscreteProcess(theInelasticProcess);
409 }
410 else if (particleName == "proton") {
411 pmanager->AddDiscreteProcess(theElasticProcess);
412 G4HadronInelasticProcess* theInelasticProcess =
413 new G4HadronInelasticProcess( "inelastic", G4Proton::Definition() );
414 theInelasticProcess->AddDataSet(theGGHNInel);
415 theInelasticProcess->RegisterMe(bertini);
416 theInelasticProcess->RegisterMe(theTheoModel);
417 pmanager->AddDiscreteProcess(theInelasticProcess);
418 }
419 else if (particleName == "anti_proton") {
420 pmanager->AddDiscreteProcess(theElasticProcess);
421 G4HadronInelasticProcess* theInelasticProcess =
422 new G4HadronInelasticProcess( "inelastic", G4AntiProton::Definition() );
423 theInelasticProcess->AddDataSet(theGGHNInel);
424 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
425 pmanager->AddDiscreteProcess(theInelasticProcess);
426
427 } else if (particleName == "neutron") {
428 // elastic scattering
429 pmanager->AddDiscreteProcess(theElasticProcess);
430
431 // inelastic scattering
432 G4HadronInelasticProcess* theInelasticProcess =
433 new G4HadronInelasticProcess( "inelastic", G4Neutron::Definition() );
434 theInelasticProcess->AddDataSet(new G4NeutronInelasticXS());
435 theInelasticProcess->RegisterMe(bertini);
436 theInelasticProcess->RegisterMe(theTheoModel);
437 pmanager->AddDiscreteProcess(theInelasticProcess);
438
439 // fission
440 G4NeutronFissionProcess* theFissionProcess = new G4NeutronFissionProcess;
441 G4LFission* theFissionModel = new G4LFission;
442 theFissionProcess->RegisterMe(theFissionModel);
443 pmanager->AddDiscreteProcess(theFissionProcess);
444
445 // capture
446 G4NeutronCaptureProcess* theCaptureProcess = new G4NeutronCaptureProcess;
447 G4NeutronRadCapture* theCaptureModel = new G4NeutronRadCapture;
448 theCaptureProcess->RegisterMe(theCaptureModel);
449 pmanager->AddDiscreteProcess(theCaptureProcess);
450
451 } else if (particleName == "anti_neutron") {
452 pmanager->AddDiscreteProcess(theElasticProcess);
453 G4HadronInelasticProcess* theInelasticProcess =
454 new G4HadronInelasticProcess( "inelastic", G4AntiNeutron::Definition() );
455 theInelasticProcess->AddDataSet(theGGHNInel);
456 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
457 pmanager->AddDiscreteProcess(theInelasticProcess);
458
459 } else if (particleName == "lambda") {
460 pmanager->AddDiscreteProcess(theElasticProcess);
461 G4HadronInelasticProcess* theInelasticProcess =
462 new G4HadronInelasticProcess( "inelastic", G4Lambda::Definition() );
463 theInelasticProcess->AddDataSet(theGGHNInel);
464 theInelasticProcess->RegisterMe(bertini);
465 theInelasticProcess->RegisterMe(theTheoModel);
466 pmanager->AddDiscreteProcess(theInelasticProcess);
467 }
468 else if (particleName == "anti_lambda") {
469 pmanager->AddDiscreteProcess(theElasticProcess);
470 G4HadronInelasticProcess* theInelasticProcess =
471 new G4HadronInelasticProcess( "inelastic", G4AntiLambda::Definition() );
472 theInelasticProcess->AddDataSet(theGGHNInel);
473 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
474 pmanager->AddDiscreteProcess(theInelasticProcess);
475 }
476 else if (particleName == "sigma+") {
477 pmanager->AddDiscreteProcess(theElasticProcess);
478 G4HadronInelasticProcess* theInelasticProcess =
479 new G4HadronInelasticProcess( "inelastic", G4SigmaPlus::Definition() );
480 theInelasticProcess->AddDataSet(theGGHNInel);
481 theInelasticProcess->RegisterMe(bertini);
482 theInelasticProcess->RegisterMe(theTheoModel);
483 pmanager->AddDiscreteProcess(theInelasticProcess);
484 }
485 else if (particleName == "sigma-") {
486 pmanager->AddDiscreteProcess(theElasticProcess);
487 G4HadronInelasticProcess* theInelasticProcess =
488 new G4HadronInelasticProcess( "inelastic", G4SigmaMinus::Definition() );
489 theInelasticProcess->AddDataSet(theGGHNInel);
490 theInelasticProcess->RegisterMe(bertini);
491 theInelasticProcess->RegisterMe(theTheoModel);
492 pmanager->AddDiscreteProcess(theInelasticProcess);
493 }
494 else if (particleName == "anti_sigma+") {
495 pmanager->AddDiscreteProcess(theElasticProcess);
496 G4HadronInelasticProcess* theInelasticProcess =
497 new G4HadronInelasticProcess( "inelastic", G4AntiSigmaPlus::Definition() );
498 theInelasticProcess->AddDataSet(theGGHNInel);
499 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
500 pmanager->AddDiscreteProcess(theInelasticProcess);
501 }
502 else if (particleName == "anti_sigma-") {
503 pmanager->AddDiscreteProcess(theElasticProcess);
504 G4HadronInelasticProcess* theInelasticProcess =
505 new G4HadronInelasticProcess( "inelastic", G4AntiSigmaMinus::Definition() );
506 theInelasticProcess->AddDataSet(theGGHNInel);
507 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
508 pmanager->AddDiscreteProcess(theInelasticProcess);
509 }
510 else if (particleName == "xi0") {
511 pmanager->AddDiscreteProcess(theElasticProcess);
512 G4HadronInelasticProcess* theInelasticProcess =
513 new G4HadronInelasticProcess( "inelastic", G4XiZero::Definition() );
514 theInelasticProcess->AddDataSet(theGGHNInel);
515 theInelasticProcess->RegisterMe(bertini);
516 theInelasticProcess->RegisterMe(theTheoModel);
517 pmanager->AddDiscreteProcess(theInelasticProcess);
518 }
519 else if (particleName == "xi-") {
520 pmanager->AddDiscreteProcess(theElasticProcess);
521 G4HadronInelasticProcess* theInelasticProcess =
522 new G4HadronInelasticProcess( "inelastic", G4XiMinus::Definition() );
523 theInelasticProcess->AddDataSet(theGGHNInel);
524 theInelasticProcess->RegisterMe(bertini);
525 theInelasticProcess->RegisterMe(theTheoModel);
526 pmanager->AddDiscreteProcess(theInelasticProcess);
527 }
528 else if (particleName == "anti_xi0") {
529 pmanager->AddDiscreteProcess(theElasticProcess);
530 G4HadronInelasticProcess* theInelasticProcess =
531 new G4HadronInelasticProcess( "inelastic", G4AntiXiZero::Definition() );
532 theInelasticProcess->AddDataSet(theGGHNInel);
533 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
534 pmanager->AddDiscreteProcess(theInelasticProcess);
535 }
536 else if (particleName == "anti_xi-") {
537 pmanager->AddDiscreteProcess(theElasticProcess);
538 G4HadronInelasticProcess* theInelasticProcess =
539 new G4HadronInelasticProcess( "inelastic", G4AntiXiMinus::Definition() );
540 theInelasticProcess->AddDataSet(theGGHNInel);
541 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
542 pmanager->AddDiscreteProcess(theInelasticProcess);
543 }
544 else if (particleName == "deuteron") {
545 pmanager->AddDiscreteProcess(theElasticProcess);
546 G4HadronInelasticProcess* theInelasticProcess =
547 new G4HadronInelasticProcess( "inelastic", G4Deuteron::Definition() );
548 theInelasticProcess->RegisterMe(binaryCascade);
549 theInelasticProcess->RegisterMe(qmd);
550 theInelasticProcess->AddDataSet(ionXS);
551 pmanager->AddDiscreteProcess(theInelasticProcess);
552 }
553 else if (particleName == "triton") {
554 pmanager->AddDiscreteProcess(theElasticProcess);
555 G4HadronInelasticProcess* theInelasticProcess =
556 new G4HadronInelasticProcess( "inelastic", G4Triton::Definition() );
557 theInelasticProcess->RegisterMe(binaryCascade);
558 theInelasticProcess->RegisterMe(qmd);
559 theInelasticProcess->AddDataSet(ionXS);
560 pmanager->AddDiscreteProcess(theInelasticProcess);
561 }
562 else if (particleName == "alpha") {
563 pmanager->AddDiscreteProcess(theElasticProcess);
564 G4HadronInelasticProcess* theInelasticProcess =
565 new G4HadronInelasticProcess( "inelastic", G4Alpha::Definition() );
566 theInelasticProcess->RegisterMe(binaryCascade);
567 theInelasticProcess->RegisterMe(qmd);
568 theInelasticProcess->AddDataSet(ionXS);
569 pmanager->AddDiscreteProcess(theInelasticProcess);
570
571 } else if (particleName == "omega-") {
572 pmanager->AddDiscreteProcess(theElasticProcess);
573 G4HadronInelasticProcess* theInelasticProcess =
574 new G4HadronInelasticProcess( "inelastic", G4OmegaMinus::Definition() );
575 theInelasticProcess->AddDataSet(theGGHNInel);
576 theInelasticProcess->RegisterMe(bertini);
577 theInelasticProcess->RegisterMe(theTheoModel);
578 pmanager->AddDiscreteProcess(theInelasticProcess);
579
580 } else if (particleName == "anti_omega-") {
581 pmanager->AddDiscreteProcess(theElasticProcess);
582 G4HadronInelasticProcess* theInelasticProcess =
583 new G4HadronInelasticProcess( "inelastic", G4AntiOmegaMinus::Definition() );
584 theInelasticProcess->AddDataSet(theGGHNInel);
585 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
586 pmanager->AddDiscreteProcess(theInelasticProcess);
587 }
588 }
589}
590
591//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
592
595
596//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
597
598#include "G4Decay.hh"
600{
601 G4Decay* theDecayProcess = new G4Decay();
602 auto particleIterator=GetParticleIterator();
603 particleIterator->reset();
604 while( (*particleIterator)() ){
605 G4ParticleDefinition* particle = particleIterator->value();
606 G4ProcessManager* pmanager = particle->GetProcessManager();
607 if (theDecayProcess->IsApplicable(*particle)) {
608 pmanager ->AddProcess(theDecayProcess);
609 pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
610 pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
611 }
612 }
613}
614
615//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
616
618{
619 if (verboseLevel >0)
620 {
621 G4cout << "B03PhysicsList::SetCuts:";
622 G4cout << "CutLength : " << defaultCutValue/mm << " (mm)" << G4endl;
623 }
624 // "G4VUserPhysicsList::SetCutsWithDefault" method sets
625 // the default cut value for all particle types
626 SetCutsWithDefault();
627}
628
629//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
630
631#include "G4ParallelWorldProcess.hh"
633
634 G4int npw = fParaWorldName.size();
635 for ( G4int i = 0; i < npw; i++){
636 G4String procName = "ParaWorldProc_"+fParaWorldName[i];
637 G4ParallelWorldProcess* theParallelWorldProcess
638 = new G4ParallelWorldProcess(procName);
639 theParallelWorldProcess->SetParallelWorld(fParaWorldName[i]);
640
641 auto particleIterator=GetParticleIterator();
642 particleIterator->reset();
643 while( (*particleIterator)() ){
644 G4ParticleDefinition* particle = particleIterator->value();
645 G4ProcessManager* pmanager = particle->GetProcessManager();
646 pmanager->AddProcess(theParallelWorldProcess);
647 if(theParallelWorldProcess->IsAtRestRequired(particle))
648 {pmanager->SetProcessOrdering(theParallelWorldProcess, idxAtRest, 9900);}
649 pmanager->SetProcessOrderingToSecond(theParallelWorldProcess, idxAlongStep);
650 pmanager->SetProcessOrdering(theParallelWorldProcess, idxPostStep, 9900);
651 }
652 }
653
654}
655
656//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
657
658#include "G4ImportanceProcess.hh"
659#include "G4IStore.hh"
661
662 G4cout << " Preparing Importance Sampling with GhostWorld "
663 << fBiasWorldName << G4endl;
664
665 G4IStore* iStore = G4IStore::GetInstance(fBiasWorldName);
666 G4GeometrySampler fGeomSampler(fBiasWorldName,"neutron");
667 fGeomSampler.SetParallel(true); // parallelworld
668 // fGeomSampler.SetWorld(iStore->GetParallelWorldVolumePointer());
669 // fGeomSampler->PrepareImportanceSampling(G4IStore::
670 // GetInstance(fBiasWorldName), 0);
671 static G4bool first = true;
672 if(first) {
673 fGeomSampler.PrepareImportanceSampling(iStore, 0);
674
675 fGeomSampler.Configure();
676 G4cout << " GeomSampler Configured!!! " << G4endl;
677 first = false;
678 }
679
680#ifdef G4MULTITHREADED
681 if(!G4Threading::IsMasterThread()) fGeomSampler.AddProcess();
682#else
683 G4cout << " Running in singlethreaded mode!!! " << G4endl;
684#endif
685
686}
687
688//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Definition of the B03PhysicsList class.
virtual void ConstructLeptHad()
virtual void ConstructGeneral()
void ConstructAllLeptons()
B03PhysicsList(G4String)
virtual void ConstructParticle()
virtual void SetCuts()
void ConstructAllBosons()
std::vector< G4String > fParaWorldName
virtual ~B03PhysicsList()
virtual void ConstructProcess()
virtual void ConstructHad()
virtual void ConstructEM()
G4String fBiasWorldName
void ConstructAllShortLiveds()

Applications | User Support | Publications | Collaboration