Geant4  10.00.p03
G4PhysicsTableHelper.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 // $Id: G4PhysicsTableHelper.cc 70369 2013-05-29 14:59:24Z gcosmo $
27 //
28 //
29 // ------------------------------------------------------------
30 // GEANT 4 class header file
31 //
32 // Class Description
33 // G4PhysicsTableHelper is a static utility class
34 // for helping proceeses to build their physics table
35 //
36 // ------------------------------------------------------------
37 // First Implementation 20 Aug. 2004 H.Kurashige
38 //
39 // ------------------------------------------------------------
40 
41 #include "G4PhysicsTableHelper.hh"
42 #include "G4ProductionCutsTable.hh"
43 
45 
47 {
48 }
49 
51 {
52 }
53 
55 {
56 }
57 
59 {
60  return *this;
61 }
62 
63 
65 {
67  size_t numberOfMCC = cutTable->GetTableSize();
68 
69  if ( physTable !=0) {
70  // compare size of physics table and number of material-cuts-couple
71  if ( physTable->size() < numberOfMCC) {
72  // enlarge physcis table
73  physTable->resize(numberOfMCC, (G4PhysicsVector*)(0));
74 #ifdef G4VERBOSE
75  if (verboseLevel>2) {
76  G4cerr << "G4PhysicsTableHelper::PreparePhysicsTable ";
77  G4cerr << "Physics Table "<< physTable ;
78  G4cerr << " is resized to " << numberOfMCC << G4endl;
79  }
80 #endif
81  } else if ( physTable->size() > numberOfMCC){
82  // ERROR: this situation should not occur
83  // size of physics table is shorter than number of material-cuts-couple
84  physTable->resize(numberOfMCC);
85 #ifdef G4VERBOSE
86  if (verboseLevel>0) {
87  G4cerr << "G4PhysicsTableHelper::PreparePhysicsTable ";
88  G4cerr << "Physics Table "<< physTable ;
89  G4cerr << " is longer than number of material-cuts-couple " << G4endl;
90  }
91 #endif
92  G4Exception( "G4PhysicsTableHelper::PreparePhysicsTable()",
93  "ProcCuts001", FatalException,
94  "Physics Table is inconsistent with material-cuts-couple");
95  }
96  } else {
97  // create PhysicsTable is given poitner is null
98  physTable = new G4PhysicsTable(numberOfMCC);
99  if (physTable!=0) {
100  physTable->resize(numberOfMCC, (G4PhysicsVector*)(0));
101  } else {
102  G4Exception( "G4PhysicsTableHelper::PreparePhysicsTable()",
103  "ProcCuts002", FatalException,
104  "Can't create Physics Table");
105  }
106  }
107 
108 #ifdef G4VERBOSE
109  if (verboseLevel>2) {
110  if ( physTable !=0) {
111  G4cerr << "Physics Table size "<< physTable->size();
112  } else {
113  G4cerr << "Physics Table does not exist ";
114  }
115  G4cerr << ": number of material-cuts-couple " << numberOfMCC << G4endl;
116  }
117 #endif
118 
119  // Reset recal-needed flag for all physics vectors
120  physTable->ResetFlagArray();
121 
122  for (size_t idx = 0; idx <numberOfMCC; idx +=1){
123  const G4MaterialCutsCouple* mcc = cutTable->GetMaterialCutsCouple(idx);
124  //check if re-calculation of the physics vector is needed
125  // MCC is not used
126  if ( !mcc->IsUsed() ) physTable->ClearFlag(idx);
127 
128  // RecalcNeeded flag of MCC is not asserted
129  if ( !mcc->IsRecalcNeeded() ) physTable->ClearFlag(idx);
130  }
131 
132  return physTable;
133 }
134 
135 
136 
138  const G4String& fileName,
139  G4bool ascii )
140 {
141  if (physTable == 0) return false;
142 
143  // retrieve physics table from the given file
144  G4PhysicsTable* tempTable = new G4PhysicsTable();
145  if (! tempTable->RetrievePhysicsTable(fileName,ascii) ){
146 #ifdef G4VERBOSE
147  if (verboseLevel>1) {
148  G4cerr << "G4PhysicsTableHelper::RetrievePhysicsTable ";
149  G4cerr << "Fail to retreive from "<< fileName << G4endl;
150  }
151 #endif
152  G4Exception( "G4ProductionCutsTable::RetrievePhysicsTable()",
153  "ProcCuts105",
154  JustWarning, "Can not retrieve physics tables from file");
155  delete tempTable;
156  return false;
157  }
158 
160  const G4MCCIndexConversionTable* converter = cutTable->GetMCCIndexConversionTable();
161 
162  // check physics table size
163  if ( tempTable->size() != converter->size()){
164 #ifdef G4VERBOSE
165  if (verboseLevel>0) {
166  G4cerr << "G4PhysicsTableHelper::RetrievePhysicsTable ";
167  G4cerr << "Size of the physics table in "<< fileName;
168  G4cerr << "( size =" << tempTable->size() << ")";
169  G4cerr << " is inconsistent with material-cut info";
170  G4cerr << "( size =" << converter->size() << ")";
171  G4cerr << G4endl;
172  }
173 #endif
174  G4Exception( "G4ProductionCutsTable::RetrievePhysicsTable()",
175  "ProcCuts106",
176  JustWarning, "Retrived file is inconsistent with current physics tables ");
177  delete tempTable;
178  return false;
179  }
180 
181  // fill the given physics table with retrived physics vectors
182  for (size_t idx=0; idx<converter->size(); idx++){
183  if (converter->IsUsed(idx)){
184  if (converter->GetIndex(idx)<0) continue;
185  size_t i = converter->GetIndex(idx);
186  G4PhysicsVector* vec = (*physTable)[i];
187  if (vec !=0 ) delete vec;
188  (*physTable)[i] = (*tempTable)[idx];
189  physTable->ClearFlag(i);
190  }
191  }
192  tempTable->clear();
193  delete tempTable;
194 
195  return true;
196 }
197 
198 
200  size_t idx,
201  G4PhysicsVector* vec)
202 {
203  if ( physTable ==0) { return; }
204 
205  if ( physTable->size() <= idx) {
206 #ifdef G4VERBOSE
207  if (verboseLevel>0) {
208  G4cerr << "G4PhysicsTableHelper::SetPhysicsVector ";
209  G4cerr << "Given index (" << idx << ") exceeds ";
210  G4cerr << "size of the physics table ";
211  G4cerr << "( size =" << physTable->size()<< ")";
212  G4cerr << G4endl;
213  }
214 #endif
215  G4Exception( "G4ProductionCutsTable::SetPhysicsVector()",
216  "ProcCuts107",
217  JustWarning, "Illegal index ");
218  return;
219  }
220 
221  // set physics vector
222  (*physTable)[idx] = vec;
223  // clear flag
224  physTable->ClearFlag(idx);
225 
226 
227 }
228 
229 
231 {
232  verboseLevel = value;
233 }
234 
236 {
237  return verboseLevel;
238 }
239 
240 
241 
static void SetPhysicsVector(G4PhysicsTable *physTable, size_t idx, G4PhysicsVector *vec)
G4bool RetrievePhysicsTable(const G4String &filename, G4bool ascii=false)
G4int GetIndex(size_t index) const
G4bool IsRecalcNeeded() const
#define G4ThreadLocal
Definition: tls.hh:52
int G4int
Definition: G4Types.hh:78
G4bool IsUsed(size_t index) const
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
static G4ThreadLocal G4int verboseLevel
G4PhysicsTableHelper & operator=(const G4PhysicsTableHelper &)
void ClearFlag(size_t i)
bool G4bool
Definition: G4Types.hh:79
const G4MCCIndexConversionTable * GetMCCIndexConversionTable() const
void resize(size_t, G4PhysicsVector *vec=(G4PhysicsVector *)(0))
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
#define G4endl
Definition: G4ios.hh:61
static G4bool RetrievePhysicsTable(G4PhysicsTable *physTable, const G4String &fileName, G4bool ascii)
G4GLOB_DLL std::ostream G4cerr
static void SetVerboseLevel(G4int value)