Geant4 Cross Reference |
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 // Please cite the following paper if you use this software 27 // Nucl.Instrum.Meth.B260:20-27, 2007 28 29 // #define MATRIX_BOUND_CHECK 30 31 #include "RunAction.hh" 32 #include "G4AnalysisManager.hh" 33 #include "G4AutoLock.hh" 34 35 namespace 36 { 37 G4Mutex aMutex = G4MUTEX_INITIALIZER; 38 } 39 40 using namespace std; 41 42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 43 44 RunAction::RunAction(DetectorConstruction* det, PrimaryGeneratorAction* pri) 45 :fDetector(det),fPrimary(pri) 46 {} 47 48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 49 50 RunAction::~RunAction() 51 {} 52 53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 54 55 void RunAction::BeginOfRunAction(const G4Run* /*aRun*/) 56 { 57 // Vector initialization 58 59 fRow=0; 60 61 fXVector = CLHEP::HepVector(32); 62 fYVector = CLHEP::HepVector(32); 63 fThetaVector = CLHEP::HepVector(32); 64 fPhiVector = CLHEP::HepVector(32); 65 66 // Histograms 67 68 // Get/create analysis manager 69 G4cout << "##### Create analysis manager " << " " << this << G4endl; 70 71 G4AnalysisManager* man = G4AnalysisManager::Instance(); 72 man->SetDefaultFileType("root"); 73 74 G4cout << "Using " << man->GetType() << " analysis manager" << G4endl; 75 76 // Open an output file 77 man->OpenFile("nanobeam"); 78 man->SetFirstHistoId(1); 79 man->SetFirstNtupleId(1); 80 81 // Create 1st ntuple (id = 1) 82 man->CreateNtuple("ntuple0", "BeamProfile"); 83 man->CreateNtupleDColumn("xIn"); 84 man->CreateNtupleDColumn("yIn"); 85 man->CreateNtupleDColumn("zIn"); 86 man->FinishNtuple(); 87 G4cout << "Ntuple-1 created" << G4endl; 88 89 // Create 2nd htuple (id = 2) 90 man->CreateNtuple("ntuple1","Grid"); 91 man->CreateNtupleDColumn("xIn"); 92 man->CreateNtupleDColumn("yIn"); 93 man->CreateNtupleDColumn("e"); 94 man->FinishNtuple(); 95 G4cout << "Ntuple-2 created" << G4endl; 96 97 // Create 3rd ntuple (id = 3) 98 man->CreateNtuple("ntuple2","Coef"); 99 man->CreateNtupleDColumn("xIn"); 100 man->CreateNtupleDColumn("yIn"); 101 man->CreateNtupleDColumn("thetaIn"); 102 man->CreateNtupleDColumn("phiIn"); 103 man->FinishNtuple(); 104 G4cout << "Ntuple-3 created" << G4endl; 105 106 return; 107 108 } 109 110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 111 112 void RunAction::EndOfRunAction(const G4Run* /*aRun*/) 113 { 114 115 if (fDetector->GetCoef()==1) 116 { 117 // 17/12/2013 - thanks to A. Dotti 118 // CLHEP Matrix inversion (as for CLHEP version 2.1.4.1) 119 // is not thread safe, we need to protect this code. 120 // Since this is not performance critical we simply lock 121 // all this part. 122 G4AutoLock l(&aMutex); 123 // 124 125 CLHEP::HepMatrix m; 126 127 // VECTOR READING 128 // VECTOR READING 129 130 m = CLHEP::HepMatrix(32,32); 131 m = fPrimary->GetMatrix(); 132 133 G4cout << G4endl; 134 G4cout << "===> NANOBEAM LINE INTRINSIC ABERRATION COEFFICIENTS (units of micrometer and mrad) :" << G4endl; 135 G4cout << G4endl; 136 137 int inv; 138 139 m.invert(inv); 140 CLHEP::HepVector tmp(32,0); 141 tmp=m*fXVector; 142 CLHEP::HepVector b; 143 b=tmp.sub(2,2); G4cout << "<x|theta>=" << b << G4endl; 144 b=tmp.sub(8,8); G4cout << "<x|theta*delta>=" << b << G4endl; 145 b=tmp.sub(10,10); G4cout << "<x|theta^3>=" << b << G4endl; 146 b=tmp.sub(12,12); G4cout << "<x|theta*phi^2>=" << b << G4endl; 147 m.invert(inv); 148 149 m.invert(inv); 150 tmp = m*fThetaVector; 151 m.invert(inv); 152 b=tmp.sub(2,2); G4cout << "<x|x>=" << b << G4endl; 153 154 m.invert(inv); 155 tmp=m*fYVector; 156 b=tmp.sub(3,3); G4cout << "<y|phi>=" << b << G4endl; 157 b=tmp.sub(9,9); G4cout << "<y|phi*delta>=" << b << G4endl; 158 b=tmp.sub(11,11); G4cout << "<y|theta^2*phi>=" << b << G4endl; 159 b=tmp.sub(13,13); G4cout << "<y|phi^3>=" << b << G4endl; 160 m.invert(inv); 161 162 m.invert(inv); 163 tmp = m*fPhiVector; 164 m.invert(inv); 165 b=tmp.sub(3,3); G4cout << "<y|y>=" << b << G4endl; 166 167 } 168 169 // Save histograms 170 171 G4AnalysisManager* man = G4AnalysisManager::Instance(); 172 man->Write(); 173 man->CloseFile(); 174 175 // Complete clean-up 176 man->Clear(); 177 } 178