Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
/*****************************************************************************
* Copyright (C) 2009 this file is part of the NPTool Project *
* *
* For the licensing terms see $NPTOOL/Licence/NPTool_Licence *
* For the list of contributors see $NPTOOL/Licence/Contributors *
*****************************************************************************/
/*****************************************************************************
* *
* Original Author : Adrien MATTA contact address: matta@ipno.in2p3.fr *
* *
* Creation Date : January 2013 *
*---------------------------------------------------------------------------*
* Decription: *
* This class deal with Beam: *
* User can enter various parameter, such as emittance or use ASCII or root *
* TH1F distribution *
* *
*---------------------------------------------------------------------------*
* Comment: *
* *
* *
*****************************************************************************/
// STL
#include <iostream>
#include <fstream>
#include <cstdlib>
#include <sstream>
#include <cmath>
#include <vector>
// NPL header
#include "NPBeam.h"
#include "NPFunction.h"
// Use CLHEP System of unit and Physical Constant
#include "CLHEP/Units/GlobalSystemOfUnits.h"
#include "CLHEP/Units/PhysicalConstants.h"
using namespace NPL;
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Beam::Beam()
{
fBeamNucleus = new Nucleus();
fEnergy = 0;
fSigmaEnergy = -1 ;
fMeanX = 0 ;
fMeanY = 0 ;
fSigmaX = -1;
fSigmaY = -0;
fMeanThetaX = 0 ;
fMeanPhiY = 0 ;
fSigmaThetaX = -1 ;
fSigmaPhiY = -1 ;
fTargetSize = 0 ;
fEffectiveTargetSize = 0 ;
fTargetThickness = 0 ;
fEffectiveTargetThickness = 0 ;
fTargetAngle = 0 ;
fTargetZ = 0 ;
// case of user given distribution
fEnergyHist = new TH1F("EnergyHist","EnergyHist",1,0,1);
fXThetaXHist = new TH2F("XThetaXHis","XThetaXHis",1,0,1,1,0,1);
fYPhiYHist = new TH2F("YPhiYHist","YPhiYHist",1,0,1,1,0,1);
}
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Beam::~Beam()
{
delete fBeamNucleus ;
delete fEnergyHist ;
delete fXThetaXHist ;
delete fYPhiYHist ;
}
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void Beam::ReadConfigurationFile(string Path){
////////General Reading needs////////
string LineBuffer;
string DataBuffer;
//////////////////////////////////////////////////////////////////////////////////////////
ifstream ReactionFile;
ReactionFile.open(Path.c_str());
bool ReadingStatus = false ;
bool check_BeamName = false ;
bool check_Energy = false ;
bool check_SigmaEnergy = false ;
bool check_MeanX = false ;
bool check_MeanY = false ;
bool check_SigmaX = false ;
bool check_SigmaY = false ;
bool check_MeanThetaX = false ;
bool check_MeanPhiY = false ;
bool check_SigmaThetaX = false ;
bool check_SigmaPhiY = false ;
bool check_EnergyProfilePath = false ;
bool check_XThetaXPath = false ;
bool check_YPhiYPath = false ;
bool check_AllEnergy = false;
bool check_AllEmittance = false;
if (!ReactionFile.is_open()) {
cout << "Error: File " << Path << " not found" << endl ;
exit(1);
}
while (!ReactionFile.eof()) {
//Pick-up next line
getline(ReactionFile, LineBuffer);
if (LineBuffer.compare(0, 4, "Beam") == 0) {
cout << "Beam Found" << endl ;
ReadingStatus = true ;
}
while(ReadingStatus){
ReactionFile >> DataBuffer;
//Search for comment Symbol %
if (DataBuffer.compare(0, 1, "%") == 0) { ReactionFile.ignore ( std::numeric_limits<std::streamsize>::max(), '\n' );}
else if (DataBuffer == "Particle=") {
check_BeamName = true ;
ReactionFile >> DataBuffer;
delete fBeamNucleus;
fBeamNucleus = new Nucleus(DataBuffer);
cout << "Beam Particle: " << fBeamNucleus->GetName() << endl;
}
else if (DataBuffer == "Energy=") {
check_Energy = true ;
ReactionFile >> DataBuffer;
fEnergy = atof(DataBuffer.c_str()) * MeV;
cout << "Beam Energy: " << fEnergy / MeV << " MeV" << endl;
}
else if (DataBuffer == "SigmaEnergy=") {
check_SigmaEnergy = true ;
ReactionFile >> DataBuffer;
fSigmaEnergy= atof(DataBuffer.c_str()) * MeV;
cout << "Beam Energy Sigma: " << fSigmaEnergy / MeV << " MeV" << endl;
}
else if (DataBuffer=="MeanX=") {
check_MeanX = true ;
ReactionFile >> DataBuffer;
fMeanX = atof(DataBuffer.c_str()) * mm;
cout << "Mean X: " << fMeanX / mm << " mm" << endl;
}
else if (DataBuffer=="MeanY=") {
check_MeanY = true ;
ReactionFile >> DataBuffer;
fMeanY = atof(DataBuffer.c_str()) * mm;
cout << "Mean Y: " << fMeanY / mm << " mm" << endl;
}
else if (DataBuffer=="SigmaX=") {
check_SigmaX = true ;
ReactionFile >> DataBuffer;
fSigmaX = atof(DataBuffer.c_str()) * mm;
cout << "Sigma X: " << fSigmaX / mm << " mm" << endl;
}
else if (DataBuffer=="SigmaY=") {
check_SigmaY = true ;
ReactionFile >> DataBuffer;
fSigmaY = atof(DataBuffer.c_str()) * mm;
cout << "Sigma Y: " << fSigmaY / mm << " mm" << endl;
}
else if (DataBuffer == "MeanThetaX=" ) {
check_MeanThetaX = true ;
ReactionFile >> DataBuffer;
fMeanThetaX = atof(DataBuffer.c_str()) * deg;
cout << "Mean Theta X: " << fMeanThetaX / deg << " deg" << endl;
}
else if (DataBuffer == "MeanPhiY=") {
check_MeanPhiY = true ;
ReactionFile >> DataBuffer;
fMeanPhiY = atof(DataBuffer.c_str()) * deg;
cout << "Mean Phi Y: " << fMeanPhiY / deg << " deg" << endl;
}
else if (DataBuffer == "SigmaThetaX=" ) {
check_SigmaThetaX = true ;
ReactionFile >> DataBuffer;
fSigmaThetaX = atof(DataBuffer.c_str()) * deg;
cout << "Sigma Theta X: " << fSigmaThetaX / deg << " deg" << endl;
}
else if (DataBuffer == "SigmaPhiY=") {
check_SigmaPhiY = true ;
ReactionFile >> DataBuffer;
fSigmaPhiY = atof(DataBuffer.c_str()) * deg;
cout << "Sigma Phi Y: " << fSigmaPhiY / deg << " deg" << endl;
}
else if (DataBuffer == "EnergyProfilePath=") {
check_EnergyProfilePath = true ;
string FileName,HistName;
ReactionFile >> FileName >> HistName;
fEnergyHist = Read1DProfile(FileName, HistName );
}
else if (DataBuffer == "XThetaXProfilePath=") {
check_XThetaXPath = true ;
string FileName,HistName;
ReactionFile >> FileName >> HistName;
fXThetaXHist = Read2DProfile(FileName, HistName );
}
else if (DataBuffer == "YPhiYProfilePath=") {
check_YPhiYPath = true ;
string FileName,HistName;
ReactionFile >> FileName >> HistName;
fYPhiYHist = Read2DProfile(FileName, HistName );
}
///////////////////////////////////////////////////
// If no Beam Token and no comment, toggle out
else
{ReadingStatus = false; cout << "WARNING : Wrong Token Sequence: Getting out " << endl ;}
///////////////////////////////////////////////////
if( ( check_MeanX && check_MeanY && check_SigmaX && check_SigmaY && check_SigmaThetaX && check_SigmaPhiY && check_MeanThetaX && check_MeanPhiY) || ( check_XThetaXPath && check_YPhiYPath ) ){
check_AllEmittance = true ;
}
if( ( check_Energy && check_SigmaEnergy ) || ( check_EnergyProfilePath ) ){
check_AllEnergy = true ;
}
// If all Token found toggle out
if( check_BeamName && check_AllEnergy && check_AllEmittance )
ReadingStatus = false ;
}
}
if( !check_BeamName || !check_AllEnergy || !check_AllEnergy )
{cout << "WARNING : Token Sequence Incomplete, Beam definition could not be Fonctionnal" << endl ;}
cout << "///////////////////////////////////////////////////" << endl << endl ;
}
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void Beam::GenerateRandomEvent(double& E, double& X, double& Y, double& Z, double& ThetaX, double& PhiY )
{
X = Y = 1000000*cm;
if(fSigmaEnergy!=-1)
E = gRandom->Gaus(fEnergy,fSigmaEnergy);
else
E = fEnergyHist->GetRandom();
if(fSigmaX!=-1){
// Shoot within the target unless target size is null (no limit)
while(sqrt(X*X+Y*Y>fEffectiveTargetSize) || fEffectiveTargetSize == 0){
NPL::RandomGaussian2D(fMeanX, fMeanThetaX, fSigmaX, fSigmaThetaX, X, ThetaX);
NPL::RandomGaussian2D(fMeanY, fMeanPhiY, fSigmaY, fSigmaPhiY, Y, PhiY);
}
}
else{
while(sqrt(X*X+Y*Y>fEffectiveTargetSize) || fEffectiveTargetSize == 0){
fXThetaXHist->GetRandom2(X,ThetaX);
fYPhiYHist->GetRandom2(Y,PhiY);
}
}
Z = fTargetZ + fEffectiveTargetThickness*(gRandom->Uniform()-0.5);
}
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void Beam::Print() const {
}