source: trunk/examples/advanced/microbeam/src/MicrobeamRunAction.cc @ 1321

Last change on this file since 1321 was 807, checked in by garnier, 16 years ago

update

  • Property svn:executable set to *
File size: 4.0 KB
Line 
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// -------------------------------------------------------------------
27// $Id: MicrobeamRunAction.cc,v 1.6 2006/06/29 16:05:37 gunter Exp $
28// -------------------------------------------------------------------
29
30#include "G4VVisManager.hh"
31#include "G4UImanager.hh"
32#include "Randomize.hh"
33
34#include "MicrobeamRunAction.hh"
35
36//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
37
38MicrobeamRunAction::MicrobeamRunAction(MicrobeamDetectorConstruction* det)
39:Detector(det)
40{   
41  saveRndm = 0; 
42}
43
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45
46MicrobeamRunAction::~MicrobeamRunAction()
47{
48  delete[] dose3DDose;
49  delete[] mapVoxels;
50}
51
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53
54void MicrobeamRunAction::BeginOfRunAction(const G4Run* /*aRun*/)
55{ 
56 
57  // save Rndm status
58  if (saveRndm > 0)
59    { 
60      CLHEP::HepRandom::showEngineStatus();
61      CLHEP::HepRandom::saveEngineStatus("beginOfRun.rndm");
62    }
63 
64  numEvent = 0;
65  nbOfHitsGas = 0;
66   
67  // ABSORBED DOSES INITIALIZATION
68  DoseN = 0;
69  DoseC = 0;
70   
71  massCytoplasm = Detector->GetMassCytoplasm();
72  massNucleus   = Detector->GetMassNucleus();
73  nbOfPixels = Detector->GetNbOfPixelsInPhantom();
74 
75  mapVoxels = new G4ThreeVector[nbOfPixels];
76  dose3DDose = new G4float[nbOfPixels];
77
78  for (G4int i=0; i<nbOfPixels; i++)
79  {
80        mapVoxels [i]=myMicrobeamPhantomConfiguration.GetVoxelThreeVector(i);
81        dose3DDose[i]=0;
82        G4ThreeVector v=mapVoxels[i];
83  }
84 
85}
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88
89void MicrobeamRunAction::EndOfRunAction(const G4Run* /*aRun*/)
90{     
91  //drawing
92  if (G4VVisManager::GetConcreteInstance()) G4UImanager::GetUIpointer()->ApplyCommand("/vis/viewer/update");
93       
94  // save Rndm status
95  if (saveRndm == 1)
96  { 
97    CLHEP::HepRandom::showEngineStatus();
98    CLHEP::HepRandom::saveEngineStatus("endOfRun.rndm");
99  }   
100 
101  FILE *myFile;
102  myFile=fopen("3DDose.txt","a");
103  for (G4int i=0; i<nbOfPixels; i++) 
104  { 
105    G4ThreeVector v;
106    v = mapVoxels[i];
107    if ( (GetNumEvent()+1) !=0) 
108      fprintf (myFile,"%f %f %f %f \n", v.x(), v.y(), v.z(), dose3DDose[i]/(GetNumEvent()+1));
109  }
110  fclose (myFile);                             
111   
112  G4cout << "-> Total number of particles detected by the gas detector : " << GetNbOfHitsGas() << G4endl; 
113  G4cout << G4endl;   
114}
Note: See TracBrowser for help on using the repository browser.