source: trunk/source/processes/hadronic/models/management/src/G4HadronicInteraction.cc @ 1348

Last change on this file since 1348 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 7.7 KB
RevLine 
[819]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//
[1315]26// $Id: G4HadronicInteraction.cc,v 1.7 2010/04/03 00:40:57 dennis Exp $
[1340]27// GEANT4 tag $Name: geant4-09-03-ref-09 $
[819]28//
[1055]29// Hadronic Interaction  base class
30// original by H.P. Wellisch
31// modified by J.L. Chuma, TRIUMF, 21-Mar-1997
32// Last modified: 04-Apr-1997
33// reimplemented 1.11.2003 JPW.
34// 23-Jan-2009 V.Ivanchenko move constructor and destructor to the body
35
[819]36#include "G4HadronicInteraction.hh"
[1055]37#include "G4HadronicInteractionRegistry.hh"
38#include "G4HadronicException.hh"
39
40G4HadronicInteraction::G4HadronicInteraction(const G4String& modelName) :
[1196]41  verboseLevel(0), theMinEnergy(0.0), theMaxEnergy(25.0*GeV), 
[1315]42  isBlocked(false), recoilEnergyThreshold(0.0), theModelName(modelName),
43  epCheckLevels(DBL_MAX, DBL_MAX)
[1055]44{ 
45  G4HadronicInteractionRegistry::Instance()->RegisterMe(this);
46}
[1315]47
48
[1055]49G4HadronicInteraction::~G4HadronicInteraction()
50{
51  G4HadronicInteractionRegistry::Instance()->RemoveMe(this);
52}
[1196]53
[1315]54
[1196]55G4double
56G4HadronicInteraction::SampleInvariantT(const G4ParticleDefinition*, 
57                                        G4double, G4int, G4int)
58{
59  return 0.0;
60}
[819]61 
[1055]62G4double G4HadronicInteraction::GetMinEnergy(
[819]63   const G4Material *aMaterial, const G4Element *anElement ) const
[1055]64{
65  size_t i;
66  if( IsBlocked(aMaterial) )return 0.*GeV;
67  if( IsBlocked(anElement) )return 0.*GeV;
68  for( i=0; i<theMinEnergyListElements.size(); ++i )
[819]69    {
[1055]70      if( anElement == theMinEnergyListElements[i].second )
71        return theMinEnergyListElements[i].first;
[819]72    }
73    for( i=0; i<theMinEnergyList.size(); ++i )
74    {
[1055]75      if( aMaterial == theMinEnergyList[i].second )
76        return theMinEnergyList[i].first;
[819]77    }
78    if(IsBlocked()) return 0.*GeV;
79    if( verboseLevel > 0 )
80      G4cout << "*** Warning from HadronicInteraction::GetMinEnergy" << G4endl
81           << "    material " << aMaterial->GetName()
82           << " not found in min energy List" << G4endl;
83    return theMinEnergy;
[1055]84}
[819]85 
[1055]86void G4HadronicInteraction::SetMinEnergy(G4double anEnergy,
[1196]87                                         const G4Element *anElement )
[1055]88{
89  if( IsBlocked(anElement) )
90    G4cout << "*** Warning from HadronicInteraction::SetMinEnergy" << G4endl
[819]91           << "    The model is not active for the Element  "
92           << anElement->GetName() << "." << G4endl;
93   
[1055]94  for( size_t i=0; i<theMinEnergyListElements.size(); ++i )
[819]95    {
96      if( anElement == theMinEnergyListElements[i].second )
97      {
98        theMinEnergyListElements[i].first = anEnergy;
99        return;
100      }
101    }
[1196]102  theMinEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
[1055]103}
[819]104 
[1055]105void G4HadronicInteraction::SetMinEnergy(G4double anEnergy,
[1196]106                                         const G4Material *aMaterial )
[1055]107{
108  if( IsBlocked(aMaterial) )
109    G4cout << "*** Warning from HadronicInteraction::SetMinEnergy" << G4endl
[819]110           << "    The model is not active for the Material "
111           << aMaterial->GetName() << "." << G4endl;
112   
[1055]113  for( size_t i=0; i<theMinEnergyList.size(); ++i )
[819]114    {
115      if( aMaterial == theMinEnergyList[i].second )
[1055]116        {
117          theMinEnergyList[i].first = anEnergy;
118          return;
119        }
[819]120    }
[1196]121  theMinEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
[1055]122}
[819]123 
[1055]124G4double G4HadronicInteraction::GetMaxEnergy(const G4Material *aMaterial, 
125                                             const G4Element *anElement ) const
126{
127  size_t i;
[1196]128  if( IsBlocked(aMaterial) )return 0.0;
129  if( IsBlocked(anElement) )return 0.0;
[1055]130  for( i=0; i<theMaxEnergyListElements.size(); ++i )
[819]131    {
[1055]132      if( anElement == theMaxEnergyListElements[i].second )
133        return theMaxEnergyListElements[i].first;
[819]134    }
[1196]135  for( i=0; i<theMaxEnergyList.size(); ++i )
[819]136    {
[1055]137      if( aMaterial == theMaxEnergyList[i].second )
138        return theMaxEnergyList[i].first;
[819]139    }
[1196]140  if(IsBlocked()) return 0.*GeV;
141  if( verboseLevel > 0 ) {
[819]142      G4cout << "*** Warning from HadronicInteraction::GetMaxEnergy" << G4endl
[1055]143             << "    material " << aMaterial->GetName()
144             << " not found in min energy List" << G4endl;
[1196]145  }
146  return theMaxEnergy;
[1055]147}
[819]148 
[1055]149void G4HadronicInteraction::SetMaxEnergy(G4double anEnergy,
[1196]150                                         const G4Element *anElement ) 
[1055]151{
[1196]152  if( IsBlocked(anElement) ) {
[1055]153    G4cout << "*** Warning from HadronicInteraction::SetMaxEnergy" << G4endl
[819]154           << "Warning: The model is not active for the Element  "
155           << anElement->GetName() << "." << G4endl;
[1196]156  }
[1055]157  for( size_t i=0; i<theMaxEnergyListElements.size(); ++i )
[819]158    {
159      if( anElement == theMaxEnergyListElements[i].second )
160      {
161        theMaxEnergyListElements[i].first = anEnergy;
162        return;
163      }
164    }
[1196]165  theMaxEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
[1055]166}
[819]167
[1055]168void G4HadronicInteraction::SetMaxEnergy(G4double anEnergy,
[1196]169                                         const G4Material *aMaterial )
[1055]170{
[1196]171  if( IsBlocked(aMaterial) ) {
[1055]172    G4cout << "*** Warning from HadronicInteraction::SetMaxEnergy" << G4endl
[819]173           << "Warning: The model is not active for the Material "
174           << aMaterial->GetName() << "." << G4endl;
[1196]175  }
[1055]176  for( size_t i=0; i<theMaxEnergyList.size(); ++i )
[819]177    {
178      if( aMaterial == theMaxEnergyList[i].second )
[1055]179        {
[1196]180          theMaxEnergyList[i].first = anEnergy;
181          return;
[1055]182        }
[819]183    }
[1196]184  theMaxEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
[1055]185}
[819]186
[1196]187void G4HadronicInteraction::DeActivateFor( const G4Material *aMaterial )
[1055]188{
189  theBlockedList.push_back(aMaterial);
190}
[819]191
[1196]192void G4HadronicInteraction::DeActivateFor( const G4Element *anElement )
[1055]193{
194  theBlockedListElements.push_back(anElement);
195}
[819]196
[1055]197G4bool G4HadronicInteraction::IsBlocked( const G4Material *aMaterial ) const
198{
199  for( size_t i=0; i<theBlockedList.size(); ++i )
[819]200    {
201      if( aMaterial == theBlockedList[i] )
[1055]202        {
203          return true;
204        }
[819]205    }
[1055]206  return false;
207}
[819]208 
[1055]209G4bool G4HadronicInteraction::IsBlocked( const G4Element *anElement ) const
210{
211  for( size_t i=0; i<theBlockedListElements.size(); ++i )
[819]212    {
213      if( anElement == theBlockedListElements[i] )
[1055]214        {
215          return true;
216        }
[819]217    }
[1055]218  return false;
219}
220
221G4HadronicInteraction::G4HadronicInteraction(const G4HadronicInteraction &right )
222{ 
223  *this = right; 
224}
225   
226const G4HadronicInteraction& 
227G4HadronicInteraction::operator=(const G4HadronicInteraction &right )
228{ 
229  G4String text = "unintended use of G4HadronicInteraction::operator=";
230  throw G4HadronicException(__FILE__, __LINE__, text); 
231  return right;
232}
[819]233 
[1055]234/* end of file */
[819]235 
Note: See TracBrowser for help on using the repository browser.