Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PhysicsModelCatalog.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// G4PhysicsModelCatalog class implementation
27//
28// Author: M.Asai (SLAC), 26 September 2013
29//
30// Revised in August 2021 by A.Ribon (CERN).
31// --------------------------------------------------------------------------
32
34#include "G4Threading.hh"
35
36G4bool G4PhysicsModelCatalog::isInitialized = false;
37std::vector< G4int >* G4PhysicsModelCatalog::theVectorOfModelIDs = nullptr;
38std::vector< G4String >* G4PhysicsModelCatalog::theVectorOfModelNames = nullptr;
39
40// --------------------------------------------------------------------------
42 if(isInitialized)
43 {
44 return;
45 }
46 if ( theVectorOfModelIDs == nullptr && theVectorOfModelNames == nullptr ) {
47 static std::vector< G4int > aVectorOfInts;
48 theVectorOfModelIDs = &aVectorOfInts;
49 static std::vector< G4String > aVectorOfStrings;
50 theVectorOfModelNames = &aVectorOfStrings;
51
52 // NOTE:
53 // The goal is that, starting from Geant4 version 11.0, all the three
54 // identifiers (modelID, index, name) remain the same regardless of the
55 // physics list, application, and version of Geant4.
56 // Therefore, after Geant4 11.0, you can only add an entry - e.g. when
57 // a new model is added, or when a pre-existing model not yet present
58 // in this catalogue is included - at the bottom of this method
59 // (rather than inserting it in the middle), and do NOT delete anything.
60 //
61 // For the modelID values, these are the considerations and choices made:
62 // - Values below 1000 are excluded because reserved to process and
63 // sub-process ID values.
64 // - Whenever resonable, modelID values should have free spaces between
65 // them, to allow for eventual, future variants - e.g. different
66 // tunings of the same model - to keep modelID values close to the
67 // original model.
68 // - modelID values are between 10'000 and 39'999, with the following
69 // subdivision in 3 categories (identified by the most significant
70 // digit):
71 // * Values between 10'000 and 19'999 are for EM models
72 // * Values between 20'000 and 29'999 are for HAD models
73 // * Values between 30'000 and 39'999 are for all the other models
74 // (i.e. included neither in EM models nor in HAD models).
75 // Note that larger values of modelID are neither more difficult to
76 // handle nor less computing performant with respect to smaller values
77 // (we remind that, for plotting, the index of the model, rather than
78 // its modelID, should be conveniently used, whereas for all the rest
79 // the modelID is recommended).
80
81 // =======================================================================
82 // ================= 1st EM MODELS : from 10'000 to 19'999 ===============
83 // =======================================================================
84
85 InsertModel( 10000, "model_EM" );
86
87 // e- production
88 InsertModel( 10010, "model_DeltaElectron" );
89 InsertModel( 10011, "model_DeltaEBelowCut" );
90 InsertModel( 10012, "model_PhotoElectron" );
91 InsertModel( 10013, "model_ComptonElectron" );
92 InsertModel( 10014, "model_TripletElectron" );
93
94 // gamma production
95 InsertModel( 10020, "model_Bremsstrahlung" );
96 InsertModel( 10021, "model_SplitBremsstrahlung" );
97 InsertModel( 10022, "model_ComptonGamma" );
98 InsertModel( 10023, "model_Annihilation" );
99 InsertModel( 10024, "model_TripletGamma" );
100 InsertModel( 10025, "model_GammaGammaEntanglement" );
101
102 // e+e- pair production
103 InsertModel( 10030, "model_EplusEminisPair" );
104
105 // atomic de-excitation
106 InsertModel( 10040, "model_Fluorescence" );
107 InsertModel( 10041, "model_gammaPIXE" );
108 InsertModel( 10050, "model_AugerElectron" );
109 InsertModel( 10051, "model_ePIXE" );
110
111 // recoil ions
112 InsertModel( 10060, "model_IonRecoil" );
113
114 // DNA models
115 InsertModel( 11000, "model_DNA" );
116 InsertModel( 11001, "model_Ritchie1994eSolvation" );
117 InsertModel( 11002, "model_Terrisol1990eSolvation" );
118 InsertModel( 11003, "model_Meesungnoen2002eSolvation" );
119 InsertModel( 11004, "model_Kreipl2009eSolvation" );
120 InsertModel( 11005, "model_MeesungnoenSolid2002eSolvation" );
121
122 // =======================================================================
123 // ============= 2nd HADRONIC MODELS : from 20'000 to 29'999 =============
124 // =======================================================================
125
126 // ----------------------------------------------------------
127 // --- Gamma- , Lepto- , Neutrino-nuclear : 20'000-20'999 ---
128 // ----------------------------------------------------------
129 // - 20'000 - 20'099 : Electromagnetic dissociation models
130 // - 20'100 - 20'199 : Gamma-nuclear models
131 // - 20'200 - 20-299 : Electron/positron-nuclear models
132 // - 20'300 - 20'399 : Muon-nuclear models
133 // - 20'400 - 20'499 : Tau-nuclear models
134 // - 20'500 - 20'999 : Neutrino models...
135 // ...
136
137 // --- EM dissociation models: 20'000 - 20'099 ---
138
139 // Class G4EMDissociation
140 InsertModel( 20000, "model_projectileEMDissociation" );
141 InsertModel( 20001, "model_targetEMDissociation" );
142
143 // --- Gamma-nuclear models: 20'100 - 20'199 ---
144
145 // Class G4LENDorBERTModel
146 InsertModel( 20100, "model_LENDorBERTModel" );
147
148 // Class G4LowEGammaNuclearModel
149 InsertModel( 20150, "model_GammaNPreco" );
150
151 // --- Charged-lepton - nuclear models: 20'200 - 20'499 ---
152
153 // Class G4ElectroVDNuclearModel
154 InsertModel( 20200, "model_G4ElectroVDNuclearModel" );
155
156 // Class G4MuonVDNuclearModel
157 InsertModel( 20300, "model_G4MuonVDNuclearModel" );
158
159 // --- Neutrino models: 20'500 - 20'999 ---
160
161 // Class G4NeutrinoElectronCcModel
162 InsertModel( 20510, "model_nu-e-inelastic" );
163
164 // Class G4NeutrinoNucleusModel
165 InsertModel( 20520, "model_neutrino-nucleus" );
166
167 // The following classes derives from G4NeutrinoNucleusModel
168 // Class G4ANuElNucleusCcModel
169 InsertModel( 20530, "model_ANuElNuclCcModel" );
170 // Class G4ANuElNucleusNcModel
171 InsertModel( 20540, "model_ANuElNuclNcModel" );
172 // Class G4ANuMuNucleusCcModel
173 InsertModel( 20550, "model_ANuMuNuclCcModel" );
174 // Class G4ANuMuNucleusNcModel
175 InsertModel( 20560, "model_ANuMuNuclNcModel" );
176 // Class G4NuElNucleusCcModel
177 InsertModel( 20570, "model_NuElNuclCcModel" );
178 // Class G4NuElNucleusNcModel
179 InsertModel( 20580, "model_NuElNuclNcModel" );
180 // Class G4NuMuNucleusCcModel
181 InsertModel( 20590, "model_NuMuNuclCcModel" );
182 // Class G4NuMuNucleusNcModel
183 InsertModel( 20600, "model_NuMuNuclNcModel" );
184
185 // ------------------------------------------------------------------------------------------
186 // --- Elastic, Charge-Exchange, Quasi-Elastic, specialized Diffraction : 21'000 - 21'999 ---
187 // ------------------------------------------------------------------------------------------
188 // - 21'000 - 21'199 : Elastic
189 // - 21'200 - 21'299 : Charge-Exchange
190 // - 21'300 - 21'499 : Quasi-Elastic
191 // - 21'500 - 21'999 : specialized Diffraction generators
192 // -----------------------------------------------------------
193
194 // --- Elastic models: 21'000 - 21'199 ---
195
196 // Class G4HadronElastic
197 InsertModel( 21000, "model_hElasticLHEP" );
198 // Class G4AntiNuclElastic
199 InsertModel( 21010, "model_AntiAElastic" );
200 // Class G4ChipsElasticModel
201 InsertModel( 21020, "model_hElasticCHIPS" );
202 // Class G4DiffuseElastic
203 InsertModel( 21030, "model_DiffuseElastic" );
204 // Class G4DiffuseElasticV2
205 InsertModel( 21040, "model_DiffuseElasticV2" );
206 // Class G4NuclNuclDiffuseElastic
207 InsertModel( 21050, "model_NNDiffuseElastic" );
208 // Class G4ElasticHadrNucleusHE
209 InsertModel( 21060, "model_hElasticGlauber" );
210 // Class G4hhElastic
211 InsertModel( 21070, "model_HadrHadrElastic" );
212 // Class G4LowEHadronElastic
213 InsertModel( 21080, "model_hLowEElastic" );
214 // Class G4LEHadronProtonElastic
215 InsertModel( 21090, "model_G4LEHadronProtonElastic" );
216 // Class G4LEnp
217 InsertModel( 21100, "model_G4LEnp" );
218 // Class G4LEpp
219 InsertModel( 21110, "model_G4LEpp" );
220 // Class G4NeutronElectronElModel
221 InsertModel( 21120, "model_n-e-elastic" );
222 // Class G4NeutrinoElectronNcModel
223 InsertModel( 21130, "model_nu-e-elastic" );
224
225 // --- Charge exchange : 21'200 - 21'299 ---
226
227 // Class: G4ChargeExchange
228 InsertModel( 21200, "model_ChargeExchange" );
229
230 // --- Quasi-Elastic : 21'300 - 21'499 ---
231
232 // Class: G4QuasiElasticChannel (which uses Chips's model G4QuasiElRatios)
233 InsertModel( 21300, "model_QuasiElastic" );
234
235 // --- Special diffraction generators : 21'500 - 21'999 ---
236
237 // Class: G4LMsdGenerator
238 InsertModel( 21500, "model_LMsdGenerator" );
239
240 // -----------------------------------------------------------------
241 // --- High energy-models (e.g. string models) : 22'000 - 22'999 ---
242 // -----------------------------------------------------------------
243 // - 22'000 - 22'099 : G4TheoFSGenerator
244 // - 22'100 - 22'199 : FTF
245 // - 22'200 - 22'299 : QGS
246 // - ...
247 // For gamma - nuclear : QGS string formation + QGS string fragmentation + Precompound
248 // For e-/e+/mu-/mu+ - nuclear : FTF string formation + Lund string fragmentation + Precompound
249 InsertModel( 22000, "model_TheoFSGenerator" );
250 // FTF string formation + Lund string fragmentation + Precompound
251 InsertModel( 22100, "model_FTFP" );
252 // FTF string formation + Lund string fragmentation + Binary Cascade
253 InsertModel( 22150, "model_FTFB" );
254 // FTF string formation + QGS string fragmentation + Precompound
255 InsertModel( 22175, "model_FTFQGSP" );
256 // QGS string formation + QGS string fragmentation + Precompound
257 InsertModel( 22200, "model_QGSP" );
258 // QGS string formation + QGS string fragmentation + Binary Cascade
259 InsertModel( 22250, "model_QGSB" );
260
261 // ----------------------------------------------------
262 // --- Intermediate energy models : 23'000 - 23'999 ---
263 // ----------------------------------------------------
264 // - 23'000 - 23'099 : BERT
265 // - 23'100 - 23'199 : BIC
266 // - 23'200 - 23'299 : INCL
267 // - 23'300 - 23'399 : QMD
268 // ...
269 // Class: G4CascadeInterface
270 InsertModel( 23000, "model_BertiniCascade" );
271
272 // The names are similar, but not identical, to those before Geant4 11.
273 // Class: G4BinaryCascade
274 InsertModel( 23100, "model_G4BinaryCascade" );
275 // Class: G4BinaryLightIonReaction
276 InsertModel( 23110, "model_G4BinaryLightIonReaction" );
277
278 // Class: G4INCLXXInterface
279 InsertModel( 23200, "model_INCLXXCascade" );
280
281 // Class: G4QMDReaction
282 InsertModel( 23300, "model_QMDModel" );
283
284 // --------------------------------------------------------------
285 // --- Pre-equilibrium/De-excitation models : 24'000 - 24'999 ---
286 // --------------------------------------------------------------
287 // - Pre-equilibrium : 24'000 - 24'099
288 // * 24'000 - 24'049 : precompound
289 // * 24'050 - 24'099 : internal BERT pre-equilibrium
290 // - de-excitation : 24'100 - 24'999
291 // * 24'100 - 24'149 : Evaporation
292 // * 24'150 - 24'199 : Photon evaporation
293 // * 24'200 - 24'299 : GEM evaporation
294 // * 24'300 - 24'349 : Fermi BreakUp
295 // * 24'350 - 24'399 : Multifragmentation
296 // * 24'400 - 24'449 : Ablation
297 // * 24'450 - 24'499 : Fission
298 // * 24'500 - 24'599 : ABLA
299 // * 24'600 - 24'699 : internal BERT de-excitation
300 // ...
301
302 // --- Pre-equilibrium: 24'000 - 24'099 ---
303
304 // Class: G4PreCompoundModel
305 InsertModel( 24000, "model_PRECO" );
306
307 // Class: G4LowEIonFragmentation
308 InsertModel( 24010, "model_LowEIonPreco" );
309
310 // Class: G4NonEquilibriumEvaporator (i.e. BERT pre-equilibrium)
311 InsertModel( 24050, "model_G4NonEquilibriumEvaporator" );
312
313 // --- De-excitation: 24'100 - 24'999 ---
314
315 // --- Evaporation: 24'100- 24'149 ---
316
317 // Class: G4EvaporationChannel
318 InsertModel( 24100, "model_G4EvaporationChannel" );
319
320 // Classes that use it: G4NeutronRadCapture and G4ExcitationHandler
321 InsertModel( 24110, "model_e-InternalConversion" );
322
323 // Class: G4UnstableFragmentBreakUp
324 InsertModel( 24120, "model_G4UnstableFragmentBreakUp" );
325
326 // --- Photon-Evaporation: 24'150 - 24'199 ---
327
328 // Class: G4PhotonEvaporation
329 InsertModel( 24150, "model_G4PhotonEvaporation" );
330
331 // Class: G4NeutronRadCapture
332 InsertModel( 24160, "model_nRadCapture" );
333
334 // --- GEM evaporation : 24'200 - 24'299 ---
335
336 // Class: G4GEMChannel
337 InsertModel( 24200, "model_G4GEMChannel" );
338
339 // Class: G4GEMChannelVI
340 InsertModel( 24210, "model_G4GEMChannelVI" );
341
342 // --- Fermi BreakUp : 24'300 - 24'349 ---
343
344 // Class: G4FermiBreakUpVI
345 InsertModel( 24300, "model_G4FermiBreakUpVI" );
346
347 // --- Multifragmentation : 24'350 - 24'399 ---
348
349 // Class: G4StatMF
350 InsertModel( 24350, "model_G4StatMF" );
351
352 // --- Ablation : 24'400 - 24'449 ---
353
354 // Class: G4WilsonAblationModel
355 InsertModel( 24400, "model_G4WilsonAblationModel" );
356
357 // --- Fission : 24'450 - 24'499 ---
358
359 // Class: G4CompetitiveFission
360 InsertModel( 24450, "model_G4CompetitiveFission" );
361
362 // --- ABLA : 24'500 - 24'599 ---
363
364 // Class: G4AblaInterface
365 InsertModel( 24500, "model_ABLA" );
366
367 // --- internal BERT de-excitation : 24'600 - 24'699 ---
368
369 // Class: G4EquilibriumEvaporator
370 InsertModel( 24600, "model_G4EquilibriumEvaporator" );
371
372 // --- Other types of de-excitation : 24'700 - 24'999 ---
373
374 // ...
375
376 // ------------------------------------------------
377 // --- Low-energy data-driven : 25'000 - 25'999 ---
378 // ------------------------------------------------
379 // - 25'000 - 25'199 : ParticleHP
380 // - 25'200 - 25'200 : LEND
381 // ...
382 // - 25'500 - 25'999 : RadioactiveDecay
383
384 // --- ParticleHP : 25'000 - 25'199 ---
385
386 // Classes: G4ParticleHPCapture , G4ParticleHPCaptureFS
387 InsertModel( 25000, "model_NeutronHPCapture" );
388
389 // Classes: G4ParticleHPElastic , G4ParticleHPElasticFS
390 InsertModel( 25010, "model_NeutronHPElastic" );
391
392 // Classes: G4ParticleHPFission , G4ParticleHPFissionFS , G4WendtFissionFragmentGenerator
393 InsertModel( 25020, "model_NeutronHPFission" );
394
395 // Inelastic: the following classes inherit from either G4ParticleHPInelasticBaseFS or
396 // G4ParticleHPInelasticCompFS, which in turn inherit from G4ParticleHPFinalState
397
398 // Class G4ParticleHPNInelasticFS
399 InsertModel( 25030, "model_G4ParticleHPNInelasticFS_F01" );
400 // Class model_G4ParticleHPNXInelasticFS
401 InsertModel( 25031, "model_G4ParticleHPNXInelasticFS_F02" );
402 // Class G4ParticleHP2NDInelasticFS
403 InsertModel( 25032, "model_G4ParticleHP2NDInelasticFS_F03" );
404 // Class G4ParticleHP2NInelasticFS
405 InsertModel( 25033, "model_G4ParticleHP2NInelasticFS_F04" );
406 // Class G4ParticleHP3NInelasticFS
407 InsertModel( 25034, "model_G4ParticleHP3NInelasticFS_F05" );
408 // Class G4ParticleHPNAInelasticFS
409 InsertModel( 25035, "model_G4ParticleHPNAInelasticFS_F06" );
410 // Class G4ParticleHPN3AInelasticFS
411 InsertModel( 25036, "model_G4ParticleHPN3AInelasticFS_F07" );
412 // Class G4ParticleHP2NAInelasticFS
413 InsertModel( 25037, "model_G4ParticleHP2NAInelasticFS_F08" );
414 // Class G4ParticleHP3NAInelasticFS
415 InsertModel( 25038, "model_G4ParticleHP3NAInelasticFS_F09" );
416 // Class G4ParticleHPNPInelasticFS
417 InsertModel( 25039, "model_G4ParticleHPNPInelasticFS_F10" );
418 // Class G4ParticleHPN2AInelasticFS
419 InsertModel( 25040, "model_G4ParticleHPN2AInelasticFS_F11" );
420 // Clas G4ParticleHP2N2AInelasticFS
421 InsertModel( 25041, "model_G4ParticleHP2N2AInelasticFS_F12" );
422 // Class G4ParticleHPNDInelasticFS
423 InsertModel( 25042, "model_G4ParticleHPNDInelasticFS_F13" );
424 // Class G4ParticleHPNTInelasticFS
425 InsertModel( 25043, "model_G4ParticleHPNTInelasticFS_F14" );
426 // Class G4ParticleHPNHe3InelasticFS
427 InsertModel( 25044, "model_G4ParticleHPNHe3InelasticFS_F15" );
428 // Class G4ParticleHPND2AInelasticFS
429 InsertModel( 25045, "model_G4ParticleHPND2AInelasticFS_F16" );
430 // Class G4ParticleHPNT2AInelasticFS
431 InsertModel( 25046, "model_G4ParticleHPNT2AInelasticFS_F17" );
432 // Class G4ParticleHP4NInelasticFS
433 InsertModel( 25047, "model_G4ParticleHP4NInelasticFS_F18" );
434 // Class G4ParticleHP2NPInelasticFS
435 InsertModel( 25048, "model_G4ParticleHP2NPInelasticFS_F19" );
436 // Class G4ParticleHP3NPInelasticFS
437 InsertModel( 25049, "model_G4ParticleHP3NPInelasticFS_F20" );
438 // Class G4ParticleHPN2PInelasticFS
439 InsertModel( 25050, "model_G4ParticleHPN2PInelasticFS_F21" );
440 // Class G4ParticleHPNPAInelasticFS
441 InsertModel( 25051, "model_G4ParticleHPNPAInelasticFS_F22" );
442 // Class G4ParticleHPPInelasticFS
443 InsertModel( 25052, "model_G4ParticleHPPInelasticFS_F23" );
444 // Class G4ParticleHPDInelasticFS
445 InsertModel( 25053, "model_G4ParticleHPDInelasticFS_F24" );
446 // Class G4ParticleHPTInelasticFS
447 InsertModel( 25054, "model_G4ParticleHPTInelasticFS_F25" );
448 // Class G4ParticleHPHe3InelasticFS
449 InsertModel( 25055, "model_G4ParticleHPHe3InelasticFS_F26" );
450 // Class G4ParticleHPAInelasticFS
451 InsertModel( 25056, "model_G4ParticleHPAInelasticFS_F27" );
452 // Class G4ParticleHP2AInelasticFS
453 InsertModel( 25057, "model_G4ParticleHP2AInelasticFS_F28" );
454 // Class G4ParticleHP3AInelasticFS
455 InsertModel( 25058, "model_G4ParticleHP3AInelasticFS_F29" );
456 // Class G4ParticleHP2PInelasticFS
457 InsertModel( 25059, "model_G4ParticleHP2PInelasticFS_F30" );
458 // Class G4ParticleHPPAInelasticFS
459 InsertModel( 25060, "model_G4ParticleHPPAInelasticFS_F31" );
460 // Class G4ParticleHPD2AInelasticFS
461 InsertModel( 25061, "model_G4ParticleHPD2AInelasticFS_F32" );
462 // Class G4ParticleHPT2AInelasticFS
463 InsertModel( 25062, "model_G4ParticleHPT2AInelasticFS_F33" );
464 // Class G4ParticleHPPDInelasticFS
465 InsertModel( 25063, "model_G4ParticleHPPDInelasticFS_F34" );
466 // Class G4ParticleHPPTInelasticFS
467 InsertModel( 25064, "model_G4ParticleHPPTInelasticFS_F35" );
468 // Class G4ParticleHPDAInelasticFS
469 InsertModel( 25065, "model_G4ParticleHPDAInelasticFS_F36" );
470
471 // --- LEND : 25'200 - 25'299 ---
472
473 // Class: G4LENDModel
474 InsertModel( 25200, "model_LENDModel" );
475
476 // The following classes inherit from G4LENDModel
477
478 // Class: G4LENDCapture
479 InsertModel( 25210, "model_LENDCapture" );
480 // Class: G4LENDElastic
481 InsertModel( 25220, "model_LENDElastic" );
482 // Class: G4LENDFission
483 InsertModel( 25230, "model_LENDFission" );
484 // Class: G4LENDInelastic
485 InsertModel( 25240, "model_LENDInelastic" );
486 // Class: G4LENDCombinedModel
487 InsertModel( 25250, "model_LENDCombinedModel" );
488 // Class: G4LENDGammaModel
489 InsertModel( 25260, "model_LENDGammaModel" );
490
491 // --- Radioactive Decay : 25'500 - 25'999 ---
492
493 // 25'510 + 10*G4RadioactiveDecayMode
494 InsertModel( 25510, "model_RDM_IT" );
495 InsertModel( 25520, "model_RDM_BetaMinus" );
496 InsertModel( 25530, "model_RDM_BetaPlus" );
497 InsertModel( 25540, "model_RDM_KshellEC" );
498 InsertModel( 25550, "model_RDM_LshellEC" );
499 InsertModel( 25560, "model_RDM_MshellEC" );
500 InsertModel( 25570, "model_RDM_NshellEC" );
501 InsertModel( 25580, "model_RDM_Alpha" );
502 InsertModel( 25590, "model_RDM_Proton" );
503 InsertModel( 25600, "model_RDM_Neutron" );
504 InsertModel( 25610, "model_RDM_SpFission" );
505 InsertModel( 25620, "model_RDM_BDProton" );
506 InsertModel( 25630, "model_RDM_BDNeutron" );
507 InsertModel( 25640, "model_RDM_Beta2Minus" );
508 InsertModel( 25650, "model_RDM_Beta2Plus" );
509 InsertModel( 25660, "model_RDM_Proton2" );
510 InsertModel( 25670, "model_RDM_Neutron2" );
511 InsertModel( 25680, "model_RDM_Triton" );
512
513 InsertModel( 25810, "model_RDM_AtomicRelaxation" );
514
515 // -------------------------------------------------------------------
516 // --- Others HAD (everything not include above) : 26'000 - 26'999 ---
517 // -------------------------------------------------------------------
518 // - 26'000 - 26'099 : Stopping
519 // - 26'100 - 26'199 : Fission
520 // - 26'200 - 26'299 : Abration
521 // - 26'300 - 26'399 : Coalescence
522 // ...
523
524 // --- Stopping : 26'000 - 26'099 ---
525
526 // Below are classes that derives from G4HadronStoppingProcess .
527 // The names are the same as before Geant4 11, except for the "model_" prefix.
528 // Classes that use it: G4HadronStoppingProcess .
529
530 // Class: G4HadronicAbsorptionBertini
531 InsertModel( 26000, "model_hBertiniCaptureAtRest_EMCascade" );
532 InsertModel( 26001, "model_hBertiniCaptureAtRest_NuclearCapture" );
533 InsertModel( 26002, "model_hBertiniCaptureAtRest_DIO" );
534 // Class: G4HadronicAbsorptionFritiof
535 InsertModel( 26010, "model_hFritiofCaptureAtRest_EMCascade" );
536 InsertModel( 26011, "model_hFritiofCaptureAtRest_NuclearCapture" );
537 InsertModel( 26012, "model_hFritiofCaptureAtRest_DIO" );
538 // Class: G4HadronicAbsorptionFritiofWithBinaryCascade
539 InsertModel( 26020, "model_hFritiofWithBinaryCascadeCaptureAtRest_EMCascade" );
540 InsertModel( 26021, "model_hFritiofWithBinaryCascadeCaptureAtRest_NuclearCapture" );
541 InsertModel( 26022, "model_hFritiofWithBinaryCascadeCaptureAtRest_DIO" );
542 // Class: G4MuonMinusCapture
543 InsertModel( 26030, "model_muMinusCaptureAtRest_EMCascade" );
544 InsertModel( 26031, "model_muMinusCaptureAtRest_NuclearCapture" );
545 InsertModel( 26032, "model_muMinusCaptureAtRest_DIO" );
546
547 // --- Fission : 26'100 - 26'199 ---
548
549 // Class G4LFission
550 InsertModel( 26100, "model_G4LFission" );
551
552 // LLNL fission (related classes: G4FissionLibrary, G4LLNLFission, G4FissLib, G4fissionEvent)
553 InsertModel( 26110, "model_G4LLNLFission" );
554
555 // --- Abration : 26'200 - 26'299 ---
556
557 // Class G4WilsonAbrasionModel
558 InsertModel( 26200, "model_G4WilsonAbrasion" );
559
560 // --- Coalescence : 26'300 - 26'399 ---
561
562 // Class G4CRCoalescence
563 InsertModel( 26300, "model_G4CRCoalescence" );
564
565 // ==========================================================================
566 // === 3rd OTHER (i.e. non-EM and non-HAD) MODELS : from 30'000 to 39'999 ===
567 // ==========================================================================
568
569 // -------------------------------
570 // --- Biasing : 30'000-30'999 ---
571 // -------------------------------
572
573 // The name is the same as before Geant4 11, except for the "model_" prefix.
574 // Classes that use it: G4BOptrForceCollision
575 InsertModel( 30010, "model_GenBiasForceCollision" );
576
577 // ----------------------------------
578 // --- Channeling : 31'000-31'999 ---
579 // ----------------------------------
580
581 // The name is the same as before Geant4 11, except for the "model_" prefix.
582 // Classes that use it: G4Channeling , G4ChannelingOptrChangeCrossSection
583 InsertModel( 31010, "model_channeling" );
584
585 // --- Others ... ---
586
587 // ======================================================================
588 // ================== 4th MODELS ADDED AFTER Geant4 11 ==================
589 // ======================================================================
590 // PLEASE ADD MODELS ONLY BELOW HERE, WITH PROPER modelID .
591 // IF YOU ARE NOT SURE, PLEASE CONTACT ONE OF THE COORDINATORS OF THE
592 // GEANT4 PHYSICS WORKING GROUPS.
593
594 // ...
595
596 SanityCheck();
597 isInitialized = true;
598
599 // The following call is commented out because it should be protected by
600 // the verbosity level, but this would imply a dependency of the global
601 // category to other categories which is not allowed.
602 // Anyhow, the call G4PhysicsModelCatalog::PrintAllInformation()
603 // can be easily placed elsewhere (in either user code, or even in other
604 // Geant4 classes).
605 //PrintAllInformation();
606 }
607}
608
609// --------------------------------------------------------------------------
610void G4PhysicsModelCatalog::SanityCheck() {
611 if ( theVectorOfModelIDs->size() != theVectorOfModelNames->size() ) {
613 ed << "theVectorOfModelIDs' size=" << theVectorOfModelIDs->size()
614 << " is NOT the same as theVectorOfModelNames's size=" << theVectorOfModelNames->size();
615 G4Exception( "G4PhysicsModelCatalog::SanityCheck()", "PhysModelCatalog001",
616 FatalException, ed, "should be the same!" );
617 } else {
618 G4bool isModelIDOutsideRange = false;
619 G4bool isModelIDRepeated = false;
620 G4bool isModelNameRepeated = false;
621 for ( int idx = 0; idx < Entries(); ++idx ) {
622 G4int modelID = (*theVectorOfModelIDs)[ idx ];
623 G4String modelName = (*theVectorOfModelNames)[ idx ];
624 if ( modelID < GetMinAllowedModelIDValue() || modelID > GetMaxAllowedModelIDValue() ) {
625 isModelIDOutsideRange = true;
626 }
627 for ( int jdx = idx + 1; jdx < Entries(); ++jdx ) {
628 if(modelID == (*theVectorOfModelIDs)[jdx])
629 {
630 isModelIDRepeated = true;
631 }
632 if(modelName == (*theVectorOfModelNames)[jdx])
633 {
634 isModelNameRepeated = true;
635 }
636 }
637 }
638 if ( isModelIDOutsideRange || isModelIDRepeated || isModelNameRepeated ) {
640 if(isModelIDOutsideRange)
641 {
642 ed << "theVectorOfModelIDs has NOT all entries between "
643 << GetMinAllowedModelIDValue() << " and "
645 }
646 if(isModelIDRepeated)
647 {
648 ed << "theVectorOfModelIDs has NOT all unique IDs !";
649 }
650 if(isModelNameRepeated)
651 {
652 ed << "theVectorOfModelNames has NOT all unique names !";
653 }
654 G4Exception( "G4PhysicsModelCatalog::SanityCheck()", "PhysModelCatalog002",
655 FatalException, ed, "cannot continue!" );
656 }
657 }
658 return;
659}
660
661// --------------------------------------------------------------------------
662
663// --------------------------------------------------------------------------
665 G4String modelName = "Undefined";
666 if ( modelID >= GetMinAllowedModelIDValue() && modelID <= GetMaxAllowedModelIDValue() ) {
667 for ( int idx = 0; idx < Entries(); ++idx ) {
668 if ( (*theVectorOfModelIDs)[ idx ] == modelID ) {
669 modelName = (*theVectorOfModelNames)[ idx ];
670 break;
671 }
672 }
673 }
674 return modelName;
675}
676
677// --------------------------------------------------------------------------
679 return ( modelIndex >= 0 && modelIndex < Entries() ) ? (*theVectorOfModelNames)[ modelIndex ] : "Undefined";
680}
681
682// --------------------------------------------------------------------------
684 return ( modelIndex >= 0 && modelIndex < Entries() ) ? (*theVectorOfModelIDs)[ modelIndex ] : -1;
685}
686
687// --------------------------------------------------------------------------
689 if(!isInitialized)
690 {
691 Initialize();
692 }
693 G4int modelID = -1;
694 for ( G4int idx = 0; idx < Entries(); ++idx ) {
695 if ( (*theVectorOfModelNames)[ idx ] == modelName ) {
696 modelID = (*theVectorOfModelIDs)[ idx ];
697 break;
698 }
699 }
700 return modelID;
701}
702
703// --------------------------------------------------------------------------
705 G4int modelIndex = -1;
706 if ( modelID >= GetMinAllowedModelIDValue() && modelID <= GetMaxAllowedModelIDValue() ) {
707 for ( G4int idx = 0; idx < Entries(); ++idx ) {
708 if ( (*theVectorOfModelIDs)[ idx ] == modelID ) {
709 modelIndex = idx;
710 break;
711 }
712 }
713 }
714 return modelIndex;
715}
716
717// --------------------------------------------------------------------------
719 G4int modelIndex = -1;
720 for ( G4int idx = 0; idx < Entries(); ++idx ) {
721 if ( (*theVectorOfModelNames)[ idx ] == modelName ) {
722 modelIndex = idx;
723 break;
724 }
725 }
726 return modelIndex;
727}
728
729// --------------------------------------------------------------------------
731 // It is enough to check the size of one of the two vectors, because they have the same size.
732 return ( theVectorOfModelIDs != nullptr ) ? G4int( theVectorOfModelIDs->size() ) : -1;
733}
734
735// --------------------------------------------------------------------------
737 G4cout << G4endl
738 << " ==================================================== " << G4endl
739 << " === G4PhysicsModelCatalog::PrintAllInformation() === " << G4endl
740 << " ==================================================== " << G4endl
741 << " SIZE (i.e. number of models in the catalog)=" << Entries() << G4endl;
742 for ( int idx = 0; idx < Entries(); ++idx ) {
743 G4int modelID = (*theVectorOfModelIDs)[ idx ];
744 G4String modelName = (*theVectorOfModelNames)[ idx ];
745 G4cout << "\t index=" << idx << "\t modelName=" << modelName
746 << "\t modelID=" << modelID << G4endl;
747 }
748 G4cout << " ==================================================== " << G4endl
749 << " ==================================================== " << G4endl
750 << " ==================================================== " << G4endl
751 << G4endl;
752}
@ 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
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static const G4String GetModelNameFromIndex(const G4int modelIndex)
static G4int GetMaxAllowedModelIDValue()
static G4int GetMinAllowedModelIDValue()
static G4int GetModelIndex(const G4int modelID)
static const G4String GetModelNameFromID(const G4int modelID)
static G4int GetModelID(const G4int modelIndex)