2 // ********************************************************************
3 // * License and Disclaimer *
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. *
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. *
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 // ********************************************************************
26 // $Id: Shielding.icc 76249 2013-11-08 11:21:29Z gcosmo $
28 //---------------------------------------------------------------------------
32 // Author: 2010 Tatsumi Koi, Gunter Folger
34 // created from FTFP_BERT
37 // 16.08.2010 H.Kurashige: Remove inclusion of G4ParticleWithCuts
38 // 26.04.2011 T.Koi: Add G4RadioactiveDecayPhysics
39 // 16.10.2012 A.Ribon: Use new default stopping
40 // 07.11.2013 T.Koi: Add IonElasticPhysics, Set proton cut to 0 to generate
41 // low energy recoils and activate production of fission
44 //----------------------------------------------------------------------------
51 #include "G4ProcessManager.hh"
52 #include "G4ProcessVector.hh"
53 #include "G4ParticleTypes.hh"
54 #include "G4ParticleTable.hh"
56 #include "G4Material.hh"
57 #include "G4MaterialTable.hh"
59 #include "G4DecayPhysics.hh"
60 #include "G4RadioactiveDecayPhysics.hh"
61 #include "G4EmProcessOptions.hh"
62 #include "G4EmStandardPhysics.hh"
63 #include "G4EmExtraPhysics.hh"
64 #include "G4IonQMDPhysics.hh"
65 #include "G4IonElasticPhysics.hh"
66 #include "G4StoppingPhysics.hh"
67 #include "G4HadronElasticPhysicsHP.hh"
68 #include "G4HadronElasticPhysicsLEND.hh"
70 #include "G4DataQuestionaire.hh"
71 #include "G4HadronPhysicsShielding.hh"
73 //template<class T> TShielding<T>::TShielding(G4int ver): T()
74 template<class T> TShielding<T>::TShielding( G4int verbose, G4String LEN_model ): T()
76 // default cut value (1.0mm)
77 // defaultCutValue = 1.0*CLHEP::mm;
78 G4DataQuestionaire it(photon, neutron, radioactive);
79 G4cout << "<<< Geant4 Physics List simulation engine: Shielding 2.1"<<G4endl;
81 this->defaultCutValue = 0.7*CLHEP::mm;
82 this->SetVerboseLevel(verbose);
85 this->RegisterPhysics( new G4EmStandardPhysics(verbose));
86 //G4EmProcessOptions emOptions;
87 //emOptions.SetFluo(true); // To activate deexcitation processes and fluorescence
88 //emOptions.SetAuger(true); // To activate Auger effect if deexcitation is activated
89 //emOptions.SetPIXE(true); // To activate Particle Induced X-Ray Emission (PIXE)
91 // Synchroton Radiation & GN Physics
92 this->RegisterPhysics( new G4EmExtraPhysics(verbose) );
95 this->RegisterPhysics( new G4DecayPhysics(verbose) );
96 //if ( rad == true ) this->RegisterPhysics( new G4RadioactiveDecayPhysics(verbose) );
97 this->RegisterPhysics( new G4RadioactiveDecayPhysics(verbose) );
99 size_t find = LEN_model.find("LEND__");
101 if ( find != G4String::npos )
103 evaluation=LEN_model;
104 evaluation.erase(0,find+6);
108 // Hadron Elastic scattering
109 if ( LEN_model == "HP" )
111 this->RegisterPhysics( new G4HadronElasticPhysicsHP(verbose) );
113 else if ( LEN_model == "LEND" )
115 this->RegisterPhysics( new G4HadronElasticPhysicsLEND(verbose,evaluation) );
116 G4DataQuestionaire itt(lend);
120 G4cout << "Shielding Physics List: Warning!" <<G4endl;
121 G4cout << "\"" << LEN_model << "\" is not valid for the low energy neutorn model." <<G4endl;
122 G4cout << "Neutron HP package will be used." <<G4endl;
123 this->RegisterPhysics( new G4HadronElasticPhysicsHP(verbose) );
127 G4HadronPhysicsShielding* hps = new G4HadronPhysicsShielding(verbose);
128 if ( LEN_model == "HP" )
132 else if ( LEN_model == "LEND" )
134 hps->UseLEND(evaluation);
138 //G4cout << "Shielding Physics List: Warning." <<G4endl;
139 //G4cout << "Name of Low Energy Neutron model " << LEN_model << " is invalid." <<G4endl;
140 //G4cout << "Will use neutron HP package." <<G4endl;
142 this->RegisterPhysics( hps );
143 //this->RegisterPhysics( new G4HadronPhysicsShielding(verbose,lend));
145 if ( LEN_model == "HP" ) {
146 //Activate prodcuton of fission fragments in neutronHP
147 char env_ff[]="G4NEUTRONHP_PRODUCE_FISSION_FRAGMENTS=1";
152 this->RegisterPhysics( new G4StoppingPhysics(verbose) );
155 this->RegisterPhysics( new G4IonQMDPhysics(verbose));
157 this->RegisterPhysics( new G4IonElasticPhysics(verbose));
159 // Neutron tracking cut --> not by default
160 // this->RegisterPhysics( new G4NeutronTrackingCut(verbose));
164 template<class T> TShielding<T>::~TShielding()
168 template<class T> void TShielding<T>::SetCuts()
170 if (this->verboseLevel >1){
171 G4cout << "Shielding::SetCuts:";
173 // " G4VUserPhysicsList::SetCutsWithDefault" method sets
174 // the default cut value for all particle types
176 this->SetCutsWithDefault();
178 //Set proton cut value to 0 for producing low energy recoil nucleus
179 this->SetCutValue(0, "proton");
181 // if (this->verboseLevel > 0)
182 // G4VUserPhysicsList::DumpCutValuesTable();