Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
exampleB03.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/exampleB03.cc
27 /// \brief Main program of the biasing/B03 example
28 //
29 //
30 // $Id: exampleB03.cc 70528 2013-05-31 16:50:36Z gcosmo $
31 //
32 //
33 // --------------------------------------------------------------
34 // GEANT 4 - exampleB03
35 //
36 // --------------------------------------------------------------
37 // Comments
38 //
39 // This example intends to show how to use both importance sampling and a
40 // customized scoring making use of the scoring framework
41 // in a parallel geometry.
42 //
43 // A simple geometry consisting of a 180 cm high concrete cylinder
44 // is constructed in the mass geometry.
45 // A geometry is constructed in the parallel geometry
46 // in order to assign importance values to slabs
47 // of width 10cm and for scoring. The parallel world volume should
48 // overlap the mass world volume and the radii of the slabs is larger
49 // than the radius of the concrete cylinder in the mass geometry.
50 // Pairs of G4GeometryCell and importance values are stored in
51 // the importance store.
52 // The scoring uses the primitive scorers via the Multi Functional Detector
53 //
54 //
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56 
57 #include <iostream>
58 
59 #ifdef G4MULTITHREADED
60 #include "G4MTRunManager.hh"
61 #else
62 #include "G4RunManager.hh"
63 #endif
64 
65 #include "G4VPhysicalVolume.hh"
66 #include "G4UImanager.hh"
67 #include "G4SystemOfUnits.hh"
68 
70 #include "B03PhysicsList.hh"
71 
73 // #include "B03PrimaryGeneratorAction.hh"
74 // #include "B03RunAction.hh"
75 
76 // construction for the parallel geometry
78 
79 // Files specific for biasing and scoring
80 //#include "G4Scorer.hh"
81 #include "G4GeometrySampler.hh"
82 #include "G4IStore.hh"
83 
84 
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87 
88 int main(int , char **)
89 {
90  G4int numberOfEvents = 100;
91 
92  G4long myseed = 345354;
93 
94 #ifdef G4MULTITHREADED
95  G4MTRunManager * runManager = new G4MTRunManager;
96  G4cout << " Number of cores: " << G4Threading::G4GetNumberOfCores() << G4endl;
97  G4cout << " but using only two! " << G4endl;
98  runManager->SetNumberOfThreads(2);
99  // runManager->SetNumberOfThreads(G4Threading::G4GetNumberOfCores());
100 #else
101  G4RunManager * runManager = new G4RunManager;
102 #endif
103 
104  G4Random::setTheSeed(myseed);
105 
106  // create the detector ---------------------------
108  runManager->SetUserInitialization(detector);
109  // ---------------------------------------------------
110 
111  // create a parallel detector
112  G4String parallelName("ParallelBiasingWorld");
114  new B03ImportanceDetectorConstruction(parallelName);
115  detector->RegisterParallelWorld(pdet);
116 
117  G4GeometrySampler pgs(pdet->GetWorldVolume(),"neutron");
118  B03PhysicsList* physlist = new B03PhysicsList;
119  physlist->AddParallelWorldName(parallelName);
120  //physlist->AddParallelWorldName(sparallelName);
121  physlist->AddBiasing(&pgs,parallelName);
122 
123  runManager->SetUserInitialization(physlist);
124 
125  // Set user action classes through Worker Initialization
126  //
128  runManager->SetUserInitialization(actions);
129 
130 // runManager->SetUserAction(new B03PrimaryGeneratorAction);
131 // // runManager->SetUserAction(new B03PrimaryGeneratorAction(ifElectron));
132 // runManager->SetUserAction(new B03RunAction);
133 
134  runManager->Initialize();
135 
136  G4VPhysicalVolume& aghostWorld = pdet->GetWorldVolumeAddress();
137  G4cout << " ghost world: " << pdet->GetName() << G4endl;
138 
139  // create an importance
140  G4IStore *aIstore = G4IStore::GetInstance(pdet->GetName());
141 
142  // create a geometry cell for the world volume replpicanumber is 0!
143  G4GeometryCell gWorldVolumeCell(aghostWorld, 0);
144  // set world volume importance to 1
145  aIstore->AddImportanceGeometryCell(1, gWorldVolumeCell);
146 
147  // set importance values and create scorers
148  G4int cell(1);
149  for (cell=1; cell<=18; cell++) {
150  G4GeometryCell gCell = pdet->GetGeometryCell(cell);
151  G4cout << " adding cell: " << cell
152  << " replica: " << gCell.GetReplicaNumber()
153  << " name: " << gCell.GetPhysicalVolume().GetName() << G4endl;
154  G4double imp = std::pow(2.0,cell-1);
155  //x aIstore.AddImportanceGeometryCell(imp, gCell);
156  aIstore->AddImportanceGeometryCell(imp, gCell.GetPhysicalVolume(), cell);
157  }
158 
159  // creating the geometry cell and add both to the store
160  // G4GeometryCell gCell = pdet->GetGeometryCell(18);
161 
162 
163  // create importance geometry cell pair for the "rest"cell
164  // with the same importance as the last concrete cell
165  G4GeometryCell gCell = pdet->GetGeometryCell(19);
166  // G4double imp = std::pow(2.0,18);
167  G4double imp = std::pow(2.0,17);
168  aIstore->AddImportanceGeometryCell(imp, gCell.GetPhysicalVolume(), 19);
169 
170  //temporary fix before runManager->BeamOn works...
171  G4UImanager* UImanager = G4UImanager::GetUIpointer();
172  G4String command1 = "/control/cout/setCoutFile fileName";
173  UImanager->ApplyCommand(command1);
174 
175  G4String command2 = "/run/beamOn "
176  + G4UIcommand::ConvertToString(numberOfEvents);;
177  UImanager->ApplyCommand(command2);
178 
179  // runManager->BeamOn(numberOfEvents);
180 
181  pgs.ClearSampling();
182 
183  delete runManager;
184 
185  return 0;
186 }
187 
188 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
virtual void SetUserInitialization(G4VUserDetectorConstruction *userInit)
void SetNumberOfThreads(G4int n)
Definition of the B03ActionInitialization class.
long G4long
Definition: G4Types.hh:80
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:357
Definition of the B03PhysicsList class.
Definition of the B03DetectorConstruction class.
int G4int
Definition: G4Types.hh:78
void RegisterParallelWorld(G4VUserParallelWorld *)
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:58
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
G4int G4GetNumberOfCores()
Definition: G4Threading.cc:102
G4int GetReplicaNumber() const
void AddImportanceGeometryCell(G4double importance, const G4GeometryCell &gCell)
Definition: G4IStore.cc:104
int main(int, char **)
Definition: exampleB03.cc:88
const G4VPhysicalVolume & GetPhysicalVolume() const
static G4IStore * GetInstance()
Definition: G4IStore.cc:211
void AddParallelWorldName(G4String &pname)
virtual void Initialize()
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Definition of the B03ImportanceDetectorConstruction class.
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:419