Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
FemtoFreezoutsAna.cxx
1/*
2 * FemtoFreezoutsAna.cxx
3 *
4 * Created on: 15-08-2014
5 * Author: Daniel Wielanek
6 * E-mail: daniel.wielanek@gmail.com
7 * Warsaw University of Technology, Faculty of Physics
8 */
9
10#include "FemtoFreezoutsAna.h"
11
12#include "Cout.h"
13#include "Cut.h"
14#include "FemtoConst.h"
15#include "FemtoFastCut.h"
16#include "FemtoFreezoutGenerator.h"
17#include "FemtoPair.h"
18#include "HistogramManager.h"
19#include "Package.h"
20#include "Parameter.h"
21
22#include <TCollection.h>
23#include <TDatabasePDG.h>
24#include <TH3.h>
25#include <TList.h>
26#include <TLorentzVector.h>
27#include <TMathBase.h>
28#include <TParticlePDG.h>
29#include <TString.h>
30#include <TVector3.h>
31
32
33namespace Hal {
34 FemtoFreezoutsAna::FemtoFreezoutsAna() :
35 TwoTrackAna(kFALSE),
36 fPdg1(211),
37 fPdg2(211),
38 fX(0),
39 fY(0),
40 fZ(0),
41 fT(0),
42 fCut(0.15),
43 fFastCut(nullptr),
44 fIgnoreSign(kFALSE),
45 fUseFakeMomenta(kFALSE),
46 fFemtoPair(nullptr),
47 fFreezoutGenerator(nullptr),
48 fHistograms1d(nullptr),
49 fHistograms3d(nullptr) {
50 fBackgroundMode = kNoBackground;
51 fMixSize = 1;
52 fKinematicsMode = EMode::kPRF;
53 for (int i = 0; i < 3; i++) {
54 fBins[i] = 100;
55 fHistoMin[i] = -100;
56 fHistoMax[i] = 100;
57 }
58 AddTags("freezouts");
59 }
60
61 FemtoFreezoutsAna::FemtoFreezoutsAna(const FemtoFreezoutsAna& ana) :
62 TwoTrackAna(ana),
63 fPdg1(ana.fPdg1),
64 fPdg2(ana.fPdg2),
65 fX(ana.fX),
66 fY(ana.fY),
67 fZ(ana.fZ),
68 fT(ana.fT),
69 fCut(ana.fCut),
70 fFastCut(nullptr),
71 fKinematicsMode(ana.fKinematicsMode),
72 fIgnoreSign(ana.fIgnoreSign),
73 fUseFakeMomenta(ana.fUseFakeMomenta),
74 fFemtoPair(nullptr),
75 fFreezoutGenerator(nullptr),
76 fHistograms1d(nullptr),
77 fHistograms3d(nullptr) {
78 if (ana.fFastCut) { fFastCut = (FemtoFastCut*) ana.fFastCut->Clone(); }
79 for (int i = 0; i < 3; i++) {
80 fBins[i] = ana.fBins[i];
81 fHistoMin[i] = ana.fHistoMin[i];
82 fHistoMax[i] = ana.fHistoMax[i];
83 }
84 if (ana.fFemtoPair) {
85 switch (ana.fKinematicsMode) {
86 case EMode::kPRF: {
87 fFemtoPair = Femto::MakePair(Femto::EKinematics::kPRF, ana.fUseFakeMomenta);
88 } break;
89 default: {
90 fFemtoPair = Femto::MakePair(Femto::EKinematics::kLCMS, ana.fUseFakeMomenta);
91 } break;
92 }
93 }
94 if (ana.fFreezoutGenerator) { fFreezoutGenerator = ana.fFreezoutGenerator->MakeCopy(); }
95 }
96
97 void FemtoFreezoutsAna::ComputePRF() {
98 Int_t bin = fFastCut->GetBin();
99 if (bin < 0) return;
100 TLorentzVector p1, p2, x1, x2;
101 p1.SetPxPyPzE(fFemtoPair->TruePx1(), fFemtoPair->TruePy1(), fFemtoPair->TruePz1(), fFemtoPair->TrueE1());
102 p2.SetPxPyPzE(fFemtoPair->TruePx2(), fFemtoPair->TruePy2(), fFemtoPair->TruePz2(), fFemtoPair->TrueE2());
103 TLorentzVector P = p1 + p2;
104 x1.SetPxPyPzE(fFemtoPair->GetX1(), fFemtoPair->GetY1(), fFemtoPair->GetZ1(), fFemtoPair->GetT1());
105 x2.SetPxPyPzE(fFemtoPair->GetX2(), fFemtoPair->GetY2(), fFemtoPair->GetZ2(), fFemtoPair->GetT2());
106 Double_t vz = P.Pz() / P.E();
107 P.Boost(0, 0, -vz); //
108 Double_t phi = P.Phi();
109 // bost in out direction
110 Double_t beta = P.Pt() / P.E();
111 p1.Boost(0, 0, -vz);
112 p1.RotateZ(-phi);
113 p1.Boost(-beta, 0, 0);
114 if (fCut < p1.P()) return;
115 // bost in out direction
116 x1.Boost(0, 0, -vz);
117 x2.Boost(0, 0, -vz);
118 // rotation
119 x1.RotateZ(-phi);
120 x2.RotateZ(-phi);
121 x1.Boost(-beta, 0, 0);
122 x2.Boost(-beta, 0, 0);
123 fX = x1.X() - x2.X();
124 fY = x1.Y() - x2.Y();
125 fZ = x1.Z() - x2.Z();
126 Double_t R = TMath::Sqrt(fX * fX + fY * fY + fZ * fZ);
127 if (fX < 0) R = -R;
128 if (fIgnoreSign) {
129 fHistograms1d->Fill(fCurrentEventCollectionID, fCurrentPairCollectionID, bin, TMath::Abs(R));
130 fHistograms3d->Fill(
131 fCurrentEventCollectionID, fCurrentPairCollectionID, bin, TMath::Abs(fX), TMath::Abs(fY), TMath::Abs(fZ));
132 } else {
133 fHistograms1d->Fill(fCurrentEventCollectionID, fCurrentPairCollectionID, bin, R);
134 fHistograms3d->Fill(fCurrentEventCollectionID, fCurrentPairCollectionID, bin, fX, fY, fZ);
135 }
136 }
137
138 void FemtoFreezoutsAna::ComputeLCMS() {
139 Int_t bin = fFastCut->GetBin();
140 if (bin < 0) return;
141 TLorentzVector p1, p2, x1, x2;
142 p1.SetPxPyPzE(fFemtoPair->TruePx1(), fFemtoPair->TruePy1(), fFemtoPair->TruePz1(), fFemtoPair->TrueE1());
143 p2.SetPxPyPzE(fFemtoPair->TruePx2(), fFemtoPair->TruePy2(), fFemtoPair->TruePz2(), fFemtoPair->TrueE2());
144 TLorentzVector P = p1 + p2;
145 TLorentzVector M = p1 - p2;
146 fT = TMath::Abs(M.M());
147 if (fT > fCut) return;
148 x1.SetPxPyPzE(fFemtoPair->GetX1(), fFemtoPair->GetY1(), fFemtoPair->GetZ1(), fFemtoPair->GetT1());
149 x2.SetPxPyPzE(fFemtoPair->GetX2(), fFemtoPair->GetY2(), fFemtoPair->GetZ2(), fFemtoPair->GetT2());
150 // boost along Z axis
151 x1.Boost(0, 0, -P.Pz() / P.E());
152 x2.Boost(0, 0, -P.Pz() / P.E());
153 P.Boost(0, 0, -P.Pz() / P.E()); //
154 // rotation
155 Double_t phi = P.Phi();
156 x1.RotateZ(-phi);
157 x2.RotateZ(-phi);
158 fX = x1.X() - x2.X();
159 fY = x1.Y() - x2.Y();
160 fZ = x1.Z() - x2.Z();
161
162 Double_t R = TMath::Sqrt(fX * fX + fY * fY + fZ * fZ);
163 if (fX < 0) R = -R;
164 if (fIgnoreSign) {
165 fHistograms1d->Fill(fCurrentEventCollectionID, fCurrentPairCollectionID, bin, TMath::Abs(R));
166 fHistograms3d->Fill(
167 fCurrentEventCollectionID, fCurrentPairCollectionID, bin, TMath::Abs(fX), TMath::Abs(fY), TMath::Abs(fZ));
168 } else {
169 fHistograms1d->Fill(fCurrentEventCollectionID, fCurrentPairCollectionID, bin, R);
170 fHistograms3d->Fill(fCurrentEventCollectionID, fCurrentPairCollectionID, bin, fX, fY, fZ);
171 }
172 }
173
174 void FemtoFreezoutsAna::ComputePRFL() {
175 Int_t bin = fFastCut->GetBin();
176 if (bin < 0) return;
177 TLorentzVector p1, p2, x1, x2;
178 p1.SetPxPyPzE(fFemtoPair->TruePx1(), fFemtoPair->TruePy1(), fFemtoPair->TruePz1(), fFemtoPair->TrueE1());
179 p2.SetPxPyPzE(fFemtoPair->TruePx2(), fFemtoPair->TruePy2(), fFemtoPair->TruePz2(), fFemtoPair->TrueE2());
180 TLorentzVector P = p1 + p2;
181 TVector3 v = P.BoostVector();
182 p1.Boost(-v);
183 if (p1.P() > fCut) return;
184 x1.SetPxPyPzE(fFemtoPair->GetX1(), fFemtoPair->GetY1(), fFemtoPair->GetZ1(), fFemtoPair->GetT1());
185 x2.SetPxPyPzE(fFemtoPair->GetX2(), fFemtoPair->GetY2(), fFemtoPair->GetZ2(), fFemtoPair->GetT2());
186 TLorentzVector X = x2 - x1;
187 X.Boost(-v);
188 // X.RotateZ(-P.Phi());
189 fX = X.X();
190 fY = X.Y();
191 fZ = X.Z();
192 fT = X.Rho();
193 Double_t R = TMath::Sqrt(fX * fX + fY * fY + fZ * fZ);
194 if (fX < 0) R = -R;
195 if (fIgnoreSign) {
196 fHistograms1d->Fill(fCurrentEventCollectionID, fCurrentPairCollectionID, bin, TMath::Abs(R));
197 fHistograms3d->Fill(
198 fCurrentEventCollectionID, fCurrentPairCollectionID, bin, TMath::Abs(fX), TMath::Abs(fY), TMath::Abs(fZ));
199 } else {
200 fHistograms1d->Fill(fCurrentEventCollectionID, fCurrentPairCollectionID, bin, R);
201 fHistograms3d->Fill(fCurrentEventCollectionID, fCurrentPairCollectionID, bin, fX, fY, fZ);
202 }
203 }
204
205 void FemtoFreezoutsAna::ComputeLCMSGamma() {
206 Int_t bin = fFastCut->GetBin();
207 if (bin < 0) return;
208 TLorentzVector p1, p2, x1, x2;
209 p1.SetPxPyPzE(fFemtoPair->TruePx1(), fFemtoPair->TruePy1(), fFemtoPair->TruePz1(), fFemtoPair->TrueE1());
210 p2.SetPxPyPzE(fFemtoPair->TruePx2(), fFemtoPair->TruePy2(), fFemtoPair->TruePz2(), fFemtoPair->TrueE2());
211 TLorentzVector P = p1 + p2;
212 TVector3 v = P.BoostVector();
213 TLorentzVector pFirst = p1;
214 pFirst.Boost(-v);
215 if (pFirst.P() > fCut) return;
216 x1.SetPxPyPzE(fFemtoPair->GetX1(), fFemtoPair->GetY1(), fFemtoPair->GetZ1(), fFemtoPair->GetT1());
217 x2.SetPxPyPzE(fFemtoPair->GetX2(), fFemtoPair->GetY2(), fFemtoPair->GetZ2(), fFemtoPair->GetT2());
218
219 Double_t vz = P.Pz() / P.E();
220 x1.Boost(0, 0, -vz);
221 x2.Boost(0, 0, -vz);
222 P.Boost(0, 0, -vz); //
223 Double_t phi = P.Phi();
224 x1.RotateZ(-phi);
225 x2.RotateZ(-phi);
226 fX = x1.X() - x2.X();
227 fY = x1.Y() - x2.Y();
228 fZ = x1.Z() - x2.Z();
229
230 p1.Boost(0, 0, -vz);
231 p1.RotateZ(-phi);
232 p2.Boost(0, 0, -vz);
233 p2.RotateZ(-phi);
234 Double_t qout = p1.X() - p2.X();
235
236
237 Double_t beta = P.Pt() / P.E();
238 x1.Boost(-beta, 0, 0);
239 x2.Boost(-beta, 0, 0);
240 p1.Boost(-beta, 0, 0);
241 Double_t kstarO = p1.X();
242 Double_t rstarO = TMath::Abs(x1.X() - x2.X());
243 fX = rstarO * kstarO / qout * 2.0;
244
245 Double_t R = TMath::Sqrt(fX * fX + fY * fY + fZ * fZ);
246 if (fX < 0) R = -R;
247 if (fIgnoreSign) {
248 fHistograms1d->Fill(fCurrentEventCollectionID, fCurrentPairCollectionID, bin, TMath::Abs(R));
249 fHistograms3d->Fill(
250 fCurrentEventCollectionID, fCurrentPairCollectionID, bin, TMath::Abs(fX), TMath::Abs(fY), TMath::Abs(fZ));
251 } else {
252 fHistograms1d->Fill(fCurrentEventCollectionID, fCurrentPairCollectionID, bin, R);
253 fHistograms3d->Fill(fCurrentEventCollectionID, fCurrentPairCollectionID, bin, fX, fY, fZ);
254 }
255 }
256
257 void FemtoFreezoutsAna::ComputeRaw() {
258 Int_t bin = fFastCut->GetBin();
259 if (bin < 0) return;
260 TLorentzVector p1, p2, x1, x2;
261 p1.SetPxPyPzE(fFemtoPair->TruePx1(), fFemtoPair->TruePy1(), fFemtoPair->TruePz1(), fFemtoPair->TrueE1());
262 p2.SetPxPyPzE(fFemtoPair->TruePx2(), fFemtoPair->TruePy2(), fFemtoPair->TruePz2(), fFemtoPair->TrueE2());
263 TLorentzVector P = p1 + p2;
264 x1.SetPxPyPzE(fFemtoPair->GetX1(), fFemtoPair->GetY1(), fFemtoPair->GetZ1(), fFemtoPair->GetT1());
265 x2.SetPxPyPzE(fFemtoPair->GetX2(), fFemtoPair->GetY2(), fFemtoPair->GetZ2(), fFemtoPair->GetT2());
266 Double_t vz = P.Pz() / P.E();
267 P.Boost(0, 0, -vz); //
268 Double_t phi = P.Phi();
269 // bost in out direction
270 Double_t beta = P.Pt() / P.E();
271 p1.Boost(0, 0, -vz);
272 p1.RotateZ(-phi);
273 p1.Boost(-beta, 0, 0);
274 if (fCut < p1.P()) return;
275 fX = x1.X() - x2.X();
276 fY = x1.Y() - x2.Y();
277 fZ = x1.Z() - x2.Z();
278 Double_t R = TMath::Sqrt(fX * fX + fY * fY + fZ * fZ);
279 if (fX < 0) R = -R;
280 if (fIgnoreSign) {
281 fHistograms1d->Fill(fCurrentEventCollectionID, fCurrentPairCollectionID, bin, TMath::Abs(R));
282 fHistograms3d->Fill(
283 fCurrentEventCollectionID, fCurrentPairCollectionID, bin, TMath::Abs(fX), TMath::Abs(fY), TMath::Abs(fZ));
284 } else {
285 fHistograms1d->Fill(fCurrentEventCollectionID, fCurrentPairCollectionID, bin, R);
286 fHistograms3d->Fill(fCurrentEventCollectionID, fCurrentPairCollectionID, bin, fX, fY, fZ);
287 }
288 }
289
290 void FemtoFreezoutsAna::ProcessFemtoPair() {
291 switch (fKinematicsMode) {
292 case EMode::kPRF: ComputePRF(); break;
293 case EMode::kLCMS: ComputeLCMS(); break;
294 case EMode::kGammaLCMS: ComputeLCMSGamma(); break;
295 case EMode::kRaw: ComputeRaw(); break;
296 default: break;
297 }
298 }
299
300 void FemtoFreezoutsAna::SetMomentumCut(Double_t cut) {
301 if (fCut <= 0) {
302 Cout::PrintInfo("KSTar cut shouldn't be <=0", EInfo::kLowWarning);
303 return;
304 }
305 fCut = cut;
306 }
307
308 void FemtoFreezoutsAna::SetOption(Option_t* opt) {
309 TString option = opt;
310 if (option.BeginsWith("background:")) {
311 Cout::PrintInfo(Form("%s backgrounds are not supported", this->ClassName()), EInfo::kDebugInfo);
312 return;
313 } else if (option.EqualTo("prf")) {
314 fKinematicsMode = EMode::kPRF;
315 } else if (option.EqualTo("lcms")) {
316 fKinematicsMode = EMode::kLCMS;
317 } else if (option.EqualTo("lcmsG")) {
318 fKinematicsMode = EMode::kGammaLCMS;
319 } else if (option.EqualTo("use_true_momenta")) {
320 fUseFakeMomenta = kFALSE;
321 } else {
322 TwoTrackAna::SetOption(opt);
323 }
324 }
325
326 void FemtoFreezoutsAna::AddCut(const Cut& cut, Option_t* opt) {
327 if (cut.GetCollectionID() != 0 && cut.GetUpdateRatio() == ECutUpdate::kEvent) {
328 Cout::PrintInfo("This class not work with more than 1 event cut collection", EInfo::kLowWarning);
329 return;
330 }
331 if (cut.GetCollectionID() >= 2 && cut.GetUpdateRatio() == ECutUpdate::kTrack) {
332 Cout::PrintInfo("This class not work with more than 1 track collections (2 "
333 "in case of non idental HBT",
334 EInfo::kLowWarning);
335 return;
336 }
337 if (cut.GetUpdateRatio() == ECutUpdate::kTwoTrack) {
338 if (cut.GetCollectionID() != 0) {
339 Cout::PrintInfo("You are trying to add two track collection with tirgger "
340 ">0, it's not supported now, try add FemtoFastCut",
341 EInfo::kLowWarning);
342 }
343 }
344 TwoTrackAna::AddCut(cut, opt);
345 }
346
347 Task::EInitFlag FemtoFreezoutsAna::Init() {
348 Task::EInitFlag prev = TwoTrackAna::Init();
349 if (prev == Task::EInitFlag::kFATAL) return prev;
350 fHistograms1d = new HistogramManager_3_1D<TH1D>();
351 fHistograms3d = new HistogramManager_3_3D<TH3D>();
352 HistogramAxisConf** axis = new HistogramAxisConf*[5];
353 if (fFastCut == NULL) { fFastCut = new FemtoFastCutVirtual(); }
354 switch (fKinematicsMode) {
355 case EMode::kPRF: {
356 fFemtoPair = Femto::MakePair(Femto::EKinematics::kPRF, fUseFakeMomenta);
357 } break;
358 default: {
359 fFemtoPair = Femto::MakePair(Femto::EKinematics::kLCMS, fUseFakeMomenta);
360 } break;
361 }
362 TParticlePDG* pid1 = TDatabasePDG::Instance()->GetParticle(fPdg1);
363 TParticlePDG* pid2 = TDatabasePDG::Instance()->GetParticle(fPdg2);
364 if (pid1) { fFemtoPair->SetMass(pid1->Mass(), pid1->Mass()); }
365 if (pid2) { fFemtoPair->SetMass(pid1->Mass(), pid2->Mass()); }
366 fFemtoPair->Init(GetTaskID());
367 fFastCut->Init(fFemtoPair);
368 TString titles[5];
369 switch (fKinematicsMode) {
370 case EMode::kPRF: {
371 titles[0] = "r^{*}_{out} [fm]";
372 titles[1] = "r^{*}_{side}[fm]";
373 titles[2] = "r^{*}_{long} [fm]";
374 titles[3] = "r^{*} [fm]";
375 } break;
376 case EMode::kLCMS: {
377 titles[0] = "r_{out} [fm]";
378 titles[1] = "r_{side}[fm]";
379 titles[2] = "r_{long} [fm]";
380 titles[3] = "r [fm]";
381 } break;
382 case EMode::kGammaLCMS: {
383 titles[0] = "r_{out #gamma} [fm]";
384 titles[1] = "r_{side}[fm]";
385 titles[2] = "r_{long} [fm]";
386 titles[3] = "r [fm]";
387 } break;
388 case EMode::kRaw: {
389 titles[0] = "#DeltaX [fm]";
390 titles[1] = "#DeltaY[fm]";
391 titles[2] = "#DeltaZ [fm]";
392 titles[3] = "#DeltaR [fm]";
393 } break;
394 }
395 titles[4] = "N_{pairs}";
396 for (int i = 0; i < 3; i++)
397 axis[i] = new HistogramAxisConf(titles[i], fBins[i], fHistoMin[i], fHistoMax[i]);
398 axis[3] = new HistogramAxisConf(titles[3], fBins[0], fHistoMin[0], fHistoMax[0]);
399 axis[4] = new HistogramAxisConf(titles[4], fBins[0], fHistoMin[0], fHistoMax[0]);
400
401 fHistograms1d->Init(fEventCollectionsNo, fTwoTrackCollectionsNo, fFastCut->GetNBins(), axis + 3, "Freezouts1d", kFALSE);
402 fHistograms3d->Init(fEventCollectionsNo, fTwoTrackCollectionsNo, fFastCut->GetNBins(), axis, "Freezouts3d", kFALSE);
403 for (int i = -0; i < 5; i++) {
404 delete axis[i];
405 }
406 delete[] axis;
407 return Task::EInitFlag::kSUCCESS;
408 }
409
410 Package* FemtoFreezoutsAna::Report() const {
411 Package* pack = TwoTrackAna::Report();
412 if (fKinematicsMode == EMode::kPRF) {
413 AddToAnaMetadata(pack, new ParameterString("FemtoKinematics", "PRF"));
414 } else if (fKinematicsMode == EMode::kLCMS) {
415 AddToAnaMetadata(pack, new ParameterString("FemtoKinematics", "LCMS"));
416 } else if (fKinematicsMode == EMode::kGammaLCMS) {
417 AddToAnaMetadata(pack, new ParameterString("FemtoKinematics", "LCMS+Gamma"));
418 }
419 if (fFreezoutGenerator) {
420 AddToAnaMetadata(pack, new ParameterString("FreezoutGenerator", "Enabled"));
421 AddToAnaMetadata(pack, fFreezoutGenerator->Report());
422 } else {
423 AddToAnaMetadata(pack, new ParameterString("FreezoutGenerator", "Disabled"));
424 }
425 if (fPdg1 != 0) {
426 AddToAnaMetadata(pack, new ParameterInt("Assumed pdg_{1}", fPdg1));
427 AddToAnaMetadata(pack, new ParameterInt("Assumed pdg_{2}", fPdg2));
428 }
429 AddToAnaMetadata(pack, new ParameterBool("IgnoreSign", fIgnoreSign));
430 TList* list = fHistograms1d->GetFlatList();
431 for (int i = 0; i < list->GetEntries(); i++) {
432 pack->AddObject(list->At(i));
433 }
434 delete list;
435 list = fHistograms3d->GetFlatList();
436 for (int i = 0; i < list->GetEntries(); i++) {
437 pack->AddObject(list->At(i));
438 }
439 delete list;
440 AddToAnaMetadata(pack, fFastCut->Report());
441 AddToAnaMetadata(pack, new ParameterDouble("freez_cut", fCut));
442 return pack;
443 }
444
445 void FemtoFreezoutsAna::ProcessPair() {
446 fFemtoPair->Build(fCurrentSignalPair);
447 PreprocessFemtoPair();
448 ProcessFemtoPair();
449 }
450
451 void FemtoFreezoutsAna::PreprocessFemtoPair() {
452 if (fFreezoutGenerator) fFreezoutGenerator->GenerateFreezoutCooordinates(fFemtoPair);
453 }
454
455 FemtoFreezoutsAna::~FemtoFreezoutsAna() {
456 delete fHistograms1d;
457 delete fHistograms3d;
458 delete fFastCut;
459 }
460
461 void FemtoFreezoutsAna::SetAxes(Int_t bins, Double_t min, Double_t max) {
462 for (int i = 0; i < 3; i++) {
463 fBins[i] = bins;
464 fHistoMin[i] = min;
465 fHistoMax[i] = max;
466 }
467 }
468
469 void FemtoFreezoutsAna::SetOutAxis(Int_t bins, Double_t min, Double_t max) {
470 const Int_t i = 0;
471 fBins[i] = bins;
472 fHistoMin[i] = min;
473 fHistoMax[i] = max;
474 }
475
476 void FemtoFreezoutsAna::SetSideAxis(Int_t bins, Double_t min, Double_t max) {
477 const Int_t i = 1;
478 fBins[i] = bins;
479 fHistoMin[i] = min;
480 fHistoMax[i] = max;
481 }
482
483 void FemtoFreezoutsAna::SetLongAxis(Int_t bins, Double_t min, Double_t max) {
484 const Int_t i = 2;
485 fBins[i] = bins;
486 fHistoMin[i] = min;
487 fHistoMax[i] = max;
488 }
489} // namespace Hal
Definition Cut.h:40
ECutUpdate GetUpdateRatio() const
Definition Cut.h:209
Int_t GetCollectionID() const
Definition Cut.h:257
Package * Report() const
Definition Package.cxx:129
void AddObject(TObject *object)
Definition Package.cxx:209