KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
KVElasticScatter.cpp
1/*
2$Id: KVElasticScatter.cpp,v 1.9 2007/04/04 10:39:17 ebonnet Exp $
3$Revision: 1.9 $
4$Date: 2007/04/04 10:39:17 $
5*/
6
7//Created by KVClassFactory on Thu Oct 19 22:31:02 2006
8//Author: John Frankland
9
10#include "KVElasticScatter.h"
11#include "KVMultiDetArray.h"
12#include "TH1F.h"
13#include "KVDetector.h"
14#include "KVTelescope.h"
15#include "KVGroup.h"
16#include "KVTarget.h"
17#include "KV2Body.h"
18#include "TObjArray.h"
19
20using namespace std;
21
23
24
25
27
28KVElasticScatter::KVElasticScatter(): fBeamDirection(0, 0, 1)
29{
30 //Default constructor
31 fDepth = 0;
32 fTheta = 0;
33 fBinE = 500;
34 fEnergy = 0;
35 fKinematics = 0;
36 fTelescope = 0;
37 fTarget = 0;
38 fIntLayer = fNDets = 0;
39 fDetector = 0;
40 fMultiLayer = kFALSE;
41 fAlignedDetectors = 0;
42 fExx = 0.;
43 fHistos = 0;
44 fDetInd = 0;
45 if (gMultiDetArray) {
46 gMultiDetArray->SetSimMode(kTRUE);
47 }
48 else {
49 Warning("KVElasticScatter", "gMultiDetArray does not refer to a valid multidetector array");
50 printf("Define it before using this class, and put it in simulation mode : gMultiDetArray->SetSimMode(kTRUE)");
51 }
52}
53
54
55
56
59
61{
62 //Destructor
63 if (fDepth)
64 delete fDepth;
65 if (fTheta)
66 delete fTheta;
67 if (fKinematics)
68 delete fKinematics;
69 if (fHistos)
70 delete fHistos;
71 if (fDetInd)
72 delete fDetInd;
73 gMultiDetArray->SetSimMode(kFALSE);
74}
75
76
77
78
81
83{
84 //Set detector parameters, target, etc. for run
85 gMultiDetArray->SetParameters(run);
86 fTarget = gMultiDetArray->GetTarget();
88 fIntLayer = 0;
90}
91
92
93
94
97
99{
100 //Set projectile Z and A
101
102 fProj.SetZandA(Z, A);
103}
104
105
106
107
110
112{
113 //Set energy of projectile in MeV
114
116 fEnergy = e;
117}
118
119
120
121
124
126{
127 //Set name of detector which will detect particle
128 fDetector = gMultiDetArray->GetDetector(det);
130 //get list of all detectors particle must pass through to get to required detector
133 //we store the association between detector type and index in list
134 if (!fDetInd)
136 else
137 fDetInd->Clear();
138 KVDetector* d;
140 Int_t i = 0;
141 while ((d = (KVDetector*) n())) {
142 //check same type is not already present
143 if (fDetInd->HasParameter(d->GetType())) {
144 //detetors with same type will be called "Type_1", "Type_2" etc.
145 TString newname;
146 int j = 1;
147 newname.Form("%s_%d", d->GetType(), j++);
148 while (fDetInd->HasParameter(newname.Data()))
149 newname.Form("%s_%d", d->GetType(), j++);
150 fDetInd->SetValue(newname.Data(), i);
151 }
152 else {
153 fDetInd->SetValue(d->GetType(), i);
154 }
155 i++;
156 }
157 //store number of aligned detectors
158 fNDets = i;
159}
160
161
162
163
169
171{
172 //For multilayer targets, use this method to choose in which
173 //layer the scattering will take place.
174 //If name="", reset any previous choice so that scattering
175 //can take place in any layer
176 if (!fTarget) {
177 cout <<
178 "<KVElasticScatter::SetTargetScatteringLayer> : No target set. Set run first."
179 << endl;
180 return;
181 }
183 if (fIntLayer)
185}
186
187
188
189
193
195{
196 //Set binning of the GetEnergy histogram
197 // Default value is 500
198 fBinE = nbins;
199}
200
201
202
206
208{
209 //Perform scattering 'N' times for current values
210 //of particle Z, A and energy, target excited state, and detector.
211
212 if (!fProj.IsDefined()) {
213 cout <<
214 "<KVElasticScatter::CalculateScattering> : Set projectile properties first"
215 << endl;
216 return;
217 }
218 if (!fEnergy) {
219 cout <<
220 "<KVElasticScatter::CalculateScattering> : Set projectile energy first"
221 << endl;
222 return;
223 }
224 if (!fDetector) {
225 cout <<
226 "<KVElasticScatter::CalculateScattering> : Set detector first" <<
227 endl;
228 return;
229 }
230 if (!fTarget) {
231 cout <<
232 "<KVElasticScatter::CalculateScattering> : No target set. Set run first."
233 << endl;
234 return;
235 }
236
237 /* -------------------------------------------------------------------------------------------------------------------------- */
238
239 //set up histograms
240
241 /* -------------------------------------------------------------------------------------------------------------------------- */
242
243 // -- histograms for debugging: depth in target and angular distribution
244 if (fDepth)
245 delete fDepth;
246 if (fTheta)
247 delete fTheta;
248 fDepth =
249 new TH1F("hDepth", "Depth (mg/cm2)", 500, 0.,
251 fTheta = new TH1F("hTheta", "Theta (deg.)", 500, 0., 0.);
252
253 /* -------------------------------------------------------------------------------------------------------------------------- */
254
255 //set up histograms for all detectors particle passes through
256 if (!fHistos) {
258 fHistos->SetOwner(); //will delete histos it stores
259 }
260 else {
261 fHistos->Clear(); //delete any previous histograms
262 }
263 KVDetector* d;
265 while ((d = (KVDetector*) n())) {
266 fHistos->
267 Add(new
268 TH1F(Form("hEloss_%s", d->GetName()), "Eloss (MeV)", fBinE, 0., 0.));
269 }
270
271 /* -------------------------------------------------------------------------------------------------------------------------- */
272
273 //set up kinematics
274 if (!fKinematics)
275 fKinematics = new KV2Body;
277 fProj.SetTheta(0);
278
279 /* -------------------------------------------------------------------------------------------------------------------------- */
280
281 //set random interaction point for scattering
282 if (fIntLayer) {
284 }
285 else {
287 }
288
289 /* -------------------------------------------------------------------------------------------------------------------------- */
290
291 //get target nucleus properties from scattering layer
293 KVMaterial* targ_mat = fTarget->GetLayer(IP);
294 KVNucleus t;
295 t.SetZ((Int_t) targ_mat->GetZ());
296 t.SetA((Int_t) targ_mat->GetMass());
298
299 /* -------------------------------------------------------------------------------------------------------------------------- */
300
301 //set excited state of target nucleus - in other words, dissipated energy for
302 //reaction due to excitation of target
304
305 /* -------------------------------------------------------------------------------------------------------------------------- */
306
307 Double_t xsec;
308 for (int i = 0; i < N; i++) {
309 //calculate slowing of incoming projectile
315 //set random direction of outgoing projectile
316
317 double th, ph;
318 th = ph = 0.;
319 fDetector->GetRandomAngles(th, ph);
320
321 //set energy of scattered nucleus
322 //WARNING: for inverse kinematics reactions, their are two energies for
323 //each angle below the maximum scattering angle.
324 //We only use the highest energy corresponding to the most forward CM angle.
325 Double_t e1, e2;
326 fKinematics->GetELab(3, th, 3, e1, e2);
327 fProj.SetEnergy(TMath::Max(e1, e2));
328 fProj.SetTheta(th);
329 fProj.SetPhi(ph);
331 fTheta->Fill(fProj.GetTheta(), xsec);
332 //slowing of outgoing projectile in target
335 //now detect particle in detector(s)
336 fAlignedDetectors->R__FOR_EACH(KVDetector, DetectParticle)(&fProj);
337 //fill histograms
338 fDepth->Fill(IP.z());
339 int j = 0;
340 n.Reset();
341 while ((d = (KVDetector*) n())) {
342 ((TH1F*)(*fHistos)[j++])->Fill(d->GetEnergy(), xsec);
343 //prepare for next round: set energy loss to zero
344 d->Clear();
345 }
347 fProj.SetTheta(0);
349 //set random interaction point for scattering
350 if (fIntLayer) {
352 }
353 else {
355 //if target is multilayer and the interaction layer is not fixed,
356 //the layer & hence the target nucleus may change
357 if (fMultiLayer) {
359 KVNucleus t;
360 t.SetZandA((Int_t) targ_mat->GetZ(), (Int_t) targ_mat->GetMass());
362 }
363 }
365 }
366}
367
368
369
370
371
377
379{
380 //Energy loss in detector of given 'type' through which scattered particle passes.
381 //Warning: if there are several detectors of the same type in the list of detectors
382 //through which the particle passes, the first one (as seen by the impinging
383 //particle) will have type "type", the second "type_1", the third "type_2", etc.
384
386}
387
388
389
390
396
398{
399 //Energy loss in any detector through which scattered particle passes.
400 //The index corresponds to the order in which detectors are crossed by the
401 //particle, beginning with 0 for the first detector, and ending with
402 //(GetNDets()-1)
403 return (TH1F*)(*fHistos)[index];
404}
405
406
407//__________________________________________________________________//
408
int Int_t
#define d(i)
#define e(i)
char Char_t
constexpr Bool_t kFALSE
double Double_t
constexpr Bool_t kTRUE
#define N
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char name[80]
char * Form(const char *fmt,...)
Relativistic binary kinematics calculator.
Definition KV2Body.h:166
void SetTarget(const KVNucleus &)
Set target for reaction.
Definition KV2Body.cpp:314
void SetEDiss(Double_t ex)
Definition KV2Body.h:246
void SetProjectile(const KVNucleus &)
Set projectile for reaction.
Definition KV2Body.cpp:339
Double_t GetXSecRuthLab(Double_t ThetaLab_Proj, Int_t OfNucleus=3) const
Definition KV2Body.cpp:1279
void SetOutgoing(const KVNucleus &proj_out)
Definition KV2Body.cpp:397
Double_t GetELab(Double_t ThetaCM, Int_t OfNucleus) const
Definition KV2Body.h:291
void CalculateKinematics()
Definition KV2Body.cpp:677
Base class for detector geometry description.
Definition KVDetector.h:160
KVGeoStrucElement * GetParentStructure(const Char_t *type, const Char_t *name="") const
void GetRandomAngles(Double_t &th, Double_t &ph, Option_t *t="isotropic")
Definition KVDetector.h:726
virtual TList * GetAlignedDetectors(UInt_t direction=1)
Calculate elastic scattering spectra in multidetector arrays.
void SetDetector(const Char_t *det)
Set name of detector which will detect particle.
Double_t fEnergy
energy of projectile
TH1F * GetEnergy()
Return pointer to energy loss histogram for chosen detector (in MeV)
Int_t fNDets
number of aligned detectors
virtual ~KVElasticScatter()
Destructor.
Int_t fIntLayer
index of interaction layer in multilayer target
Double_t fExx
excited state of target nucleus
TH1F * fDepth
depth of scattering point in target
TList * fAlignedDetectors
all aligned detectors
TVector3 fBeamDirection
beam direction vector
KVDetector * fDetector
detector where particle will be detected
KVNameValueList * fDetInd
detector type-index association
KVTelescope * fTelescope
telescope where particle will be detected
void CalculateScattering(Int_t N)
void SetEbinning(Int_t nbins=500)
void SetProjectile(Int_t Z, Int_t A)
Set projectile Z and A.
void SetRun(Int_t run)
Set detector parameters, target, etc. for run.
Bool_t fMultiLayer
kTRUE for multilayer target
void SetTargetScatteringLayer(const Char_t *name)
Int_t fBinE
Number of bins of the Energy histogram.
KV2Body * fKinematics
kinematics calculation
TH1F * fTheta
angle of scattered particle
void SetEnergy(Double_t E)
Set energy of projectile in MeV.
KVTarget * fTarget
target for current run
TObjArray * fHistos
energy loss histograms for all hit detectors
KVNucleus fProj
scattered nucleus
KVDetector * GetDetector(const Char_t *name) const
Return detector in this structure with given name.
@ kForwards
Definition KVGroup.h:33
Description of physical materials used to construct detectors & targets; interface to range tables.
Definition KVMaterial.h:94
Double_t GetZ() const
Double_t GetMass() const
virtual void SetParameters(UInt_t n, Bool_t physics_parameters_only=kFALSE)
virtual void SetSimMode(Bool_t on=kTRUE)
KVTarget * GetTarget()
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
Int_t GetIntValue(const Char_t *name) const
void SetValue(const Char_t *name, value_type value)
virtual void Clear(Option_t *opt="")
Bool_t HasParameter(const Char_t *name) const
Description of properties and kinematics of atomic nuclei.
Definition KVNucleus.h:126
Bool_t IsDefined() const
Definition KVNucleus.h:202
void SetA(Int_t a)
void SetZ(Int_t z, Char_t mt=-1)
void SetZandA(Int_t z, Int_t a)
Set atomic number and mass number.
void SetTheta(Double_t theta)
Definition KVParticle.h:693
void SetPhi(Double_t phi)
Definition KVParticle.h:697
Double_t GetTheta() const
Definition KVParticle.h:680
KVNameValueList * GetParameters() const
Definition KVParticle.h:815
void SetEnergy(Double_t e)
Definition KVParticle.h:599
void SetIncoming(Bool_t r=kTRUE)
Definition KVTarget.h:217
virtual void DetectParticle(KVNucleus *, TVector3 *norm=0)
Definition KVTarget.cpp:485
Int_t NumberOfLayers() const
Definition KVTarget.h:166
TVector3 & GetInteractionPoint(KVParticle *part=0)
Definition KVTarget.cpp:767
void SetRandomized(Bool_t r=kTRUE)
Definition KVTarget.h:250
KVMaterial * GetLayer(TVector3 &depth)
Definition KVTarget.cpp:272
Double_t GetTotalEffectiveThickness(KVParticle *part=0)
Definition KVTarget.cpp:412
void SetInteractionLayer(Int_t ilayer, TVector3 &dir)
Definition KVTarget.cpp:809
void SetOutgoing(Bool_t r=kTRUE)
Definition KVTarget.h:237
Int_t GetLayerIndex(TVector3 &depth)
Definition KVTarget.cpp:292
Associates two detectors placed one behind the other.
Definition KVTelescope.h:36
virtual void SetOwner(Bool_t enable=kTRUE)
virtual Int_t GetSize() const
virtual Int_t Fill(const char *name, Double_t w)
void Clear(Option_t *option="") override
const char * Data() const
void Form(const char *fmt,...)
Double_t z() const
const Int_t n
void Warning(const char *location, const char *fmt,...)
Add
Double_t Abs(Double_t d)
Double_t Max(Double_t a, Double_t b)
ClassImp(TPyArg)