R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandClusterMon.cxx
Go to the documentation of this file.
1/******************************************************************************
2 * Copyright (C) 2019 GSI Helmholtzzentrum für Schwerionenforschung GmbH *
3 * Copyright (C) 2019-2025 Members of R3B Collaboration *
4 * *
5 * This software is distributed under the terms of the *
6 * GNU General Public Licence (GPL) version 3, *
7 * copied verbatim in the file "LICENSE". *
8 * *
9 * In applying this license GSI does not waive the privileges and immunities *
10 * granted to it by virtue of its status as an Intergovernmental Organization *
11 * or submit itself to any jurisdiction. *
12 ******************************************************************************/
13
15#include "ElasticScattering.h"
16#include "FairLogger.h"
17#include "FairRootManager.h"
18#include "R3BNeulandCluster.h"
19#include "TClonesArray.h"
20#include "TDirectory.h"
21#include "TH1D.h"
22#include "TH2D.h"
23#include "TH3D.h"
24#include <TFile.h>
25#include <algorithm>
26#include <iostream>
27#include <numeric>
28#include <utility>
29
30namespace
31{
32 constexpr double rad2deg = 180. / 3.141592653589793238463;
33
34 inline auto GetTheta(const R3BNeulandCluster& cluster) -> double
35 {
36 const auto direction = cluster.GetLastHit().GetPosition() - cluster.GetFirstHit().GetPosition();
37 const auto x = std::acos(direction.Y() / direction.Mag()) * rad2deg;
38 // Not sure, but Kondos Theta is -90:90
39 return x - 90.;
40 }
41} // namespace
42
43R3BNeulandClusterMon::R3BNeulandClusterMon(TString input, TString output, const Option_t* option)
44 : FairTask("R3B NeuLAND NeulandCluster Monitor")
45 , fNeulandClusters(input)
46 , fOutput(std::move(output))
47 , fBeta(0.793) // 600 Mev?
48{
49 LOG(info) << "Using R3B NeuLAND NeulandCluster Monitor";
50
51 TString opt = option;
52 opt.ToUpper();
53
54 if (opt.Contains("3DTRACK"))
55 {
56 fIs3DTrackEnabled = true;
57 LOG(info) << "... with 3D track visualization";
58 }
59 else
60 {
61 fIs3DTrackEnabled = false;
62 }
63}
64
66{
67 fNeulandClusters.init();
68
69 FairRootManager* ioman = FairRootManager::Instance();
70 if (!ioman)
71 {
72 LOG(fatal) << "R3BNeulandClusterMon::Init: No FairRootManager";
73 return kFATAL;
74 }
75
76 if (fIs3DTrackEnabled)
77 {
78 // XYZ -> ZXY (side view)
79 fh3 = new TH3D("hClusters", "hClusters", 60, 1400, 1700, 50, -125, 125, 50, -125, 125);
80 fh3->SetTitle("NeuLAND Clusters");
81 fh3->GetXaxis()->SetTitle("Z");
82 fh3->GetYaxis()->SetTitle("X");
83 fh3->GetZaxis()->SetTitle("Y");
84
85 ioman->Register(fOutput, "Cluster TH3Ds in NeuLAND", fh3, kTRUE);
86 }
87
88 TH1::AddDirectory(kFALSE);
89
90 fhClusters = new TH1D("nClusters", "Number of clusters in one event", 100, 0, 100);
91 fhClusterSize = new TH1D("ClusterSize", "Number of Digis for each Cluster", 100, 0, 100);
92 fhClusterEnergy = new TH1D("ClusterEnergy", "Energy for each Cluster", 2000, 0., 2000.);
93 fhClusterEnergyVSSize =
94 new TH2D("ClusterEnergyVSSize", "Energy for each Cluster vs Cluster Size", 2000, 0., 2000., 100, 0., 100.);
95 fhClusterTime = new TH1D("ClusterTime", "Time for each Cluster", 5000, 0., 500.);
96 fhClusterEToF = new TH1D("ClusterEToF", "Cluster EToF", 2000, 0, 2000);
97 fhClusterRValue = new TH1D("ClusterRValue", "Log Cluster R Value", 20000, -10, 1);
98
99 fhClusterSizeVSEToF =
100 new TH2D("ClusterSizeVSEToF", "Number of Digis for each Cluster vs EToF", 100, 0, 100, 2000, 0, 2000);
101 fhClusterEnergyVSEToF =
102 new TH2D("ClusterEnergyVSEToF", "Energy for each Cluster vs EToF", 2000, 0., 2000., 2000, 0, 2000);
103 fhClusterEnergyVSSizeVSEToF = new TH3D("ClusterEnergyVSSizeVSEToF",
104 "Energy for each Cluster vs Cluster Size vs EToF",
105 1000,
106 0.,
107 1000.,
108 100,
109 0.,
110 100.,
111 1200,
112 0,
113 1200.);
114
115 fhENFromScatterVSEToF = new TH2D("ENFromScatterVSEToF", "ENFromScatterVSEToF", 1000, 0, 1000, 1000, 0, 1000);
116
117 fhClusterNumberVSEnergy =
118 new TH2D("ClusterNumberEtot", "Number of Clusters vs. Total Energy", 200, 0, 2000, 50, 0, 50);
119 fhClusterNumberVSEnergy->GetXaxis()->SetTitle("Total energy [MeV]");
120 fhClusterNumberVSEnergy->GetYaxis()->SetTitle("Number of Clusters");
121
122 fhClusterEToFVSEnergy =
123 new TH2D("ClusterEToFVSEnergy", "Cluster E_{ToF} vs. Cluster Energy", 100, 0, 1000, 100, 0, 1000);
124 fhClusterEToFVSEnergy->GetXaxis()->SetTitle("Cluster E_{ToF} [MeV]");
125 fhClusterEToFVSEnergy->GetYaxis()->SetTitle("Cluster E [MeV]");
126
127 fhClusterEToFVSTime = new TH2D("ClusterEToFVSTime", "Cluster E_{ToF} vs. Cluster Time", 100, 0, 1000, 500, 0, 500);
128 fhClusterEToFVSTime->GetXaxis()->SetTitle("Cluster E_{ToF} [MeV]");
129 fhClusterEToFVSTime->GetYaxis()->SetTitle("Cluster t [ns]");
130
131 fhClusterEVSTime = new TH2D("ClusterEVSTime", "Cluster E vs. Cluster Time", 250, 0, 250, 500, 0, 250);
132 fhClusterEVSTime->GetXaxis()->SetTitle("Cluster E [MeV]");
133 fhClusterEVSTime->GetYaxis()->SetTitle("Cluster t [ns]");
134
135 fhClusterLastMinusFirstDigiMagVSEnergy = new TH2D(
136 "ClusterLastMinusFirstDigiMagVSEnergy", "ClusterLastMinusFirstDigiMagVSEnergy", 101, 0, 100, 60, 0, 600);
137
138 fhClusterForemostMinusCentroidVSEnergy =
139 new TH2D("ClusterForemostMinusCentroidVSEnergy",
140 "Distance between Foremost Digi and Energy Centroid vs Cluster Energy",
141 101,
142 0,
143 100,
144 60,
145 0,
146 600);
147 fhClusterForemostMinusCentroidVSEnergy->GetXaxis()->SetTitle("|#vec{r_{z_{min}}}-#vec{r_{EC}}| [cm]");
148 fhClusterForemostMinusCentroidVSEnergy->GetYaxis()->SetTitle("Cluster Energy [MeV]");
149
150 fhClusterForemostMinusMaxEnergyDigiPosVSEnergy =
151 new TH2D("ClusterForemostMinusMaxEnergyDigiPosVSEnergy",
152 "Distance between Foremost Digi and Max Energy Digi vs Cluster Energy",
153 101,
154 0,
155 100,
156 60,
157 0,
158 600);
159 fhClusterForemostMinusMaxEnergyDigiPosVSEnergy->GetXaxis()->SetTitle("|#vec{r_{z_{min}}}-#vec{r_{E_{max}}}| [cm]");
160 fhClusterForemostMinusMaxEnergyDigiPosVSEnergy->GetYaxis()->SetTitle("Cluster Energy [MeV]");
161
162 fhClusterCentroidMinusFirstDigiPosVSEnergy =
163 new TH2D("ClusterCentroidMinusFirstDigiPosVSEnergy",
164 "Distance between Cluster Energy Centroid and First Hit vs Cluster Energy",
165 101,
166 0,
167 100,
168 60,
169 0,
170 600);
171 fhClusterCentroidMinusFirstDigiPosVSEnergy->GetXaxis()->SetTitle("|#vec{r_{EC}}-#vec{r_{t_{min}}}| [cm]");
172 fhClusterCentroidMinusFirstDigiPosVSEnergy->GetYaxis()->SetTitle("Cluster Energy [MeV]");
173
174 fhClusterMaxEnergyDigiMinusFirstDigiPosVSEnergy =
175 new TH2D("ClusterMaxEnergyDigiMinusFirstDigiPosVSEnergy",
176 "Distance between Max Energy Digi and First Digi vs Cluster Energy",
177 101,
178 0,
179 100,
180 60,
181 0,
182 600);
183 fhClusterMaxEnergyDigiMinusFirstDigiPosVSEnergy->GetXaxis()->SetTitle("|#vec{r_{E_{max}}}-#vec{r_{t_{min}}}| [cm]");
184 fhClusterMaxEnergyDigiMinusFirstDigiPosVSEnergy->GetYaxis()->SetTitle("Cluster Energy [MeV]");
185
186 fhClusterMaxEnergyDigiMinusCentroidVSEnergy =
187 new TH2D("ClusterMaxEnergyDigiMinusCentroidVSEnergy",
188 "Distance between Max Energy Digi and Energy Centroid vs Cluster Energy",
189 101,
190 0,
191 100,
192 60,
193 0,
194 600);
195 fhClusterMaxEnergyDigiMinusCentroidVSEnergy->GetXaxis()->SetTitle("|#vec{r_{E_{max}}}-#vec{r_{EC}}| [cm]");
196 fhClusterMaxEnergyDigiMinusCentroidVSEnergy->GetYaxis()->SetTitle("Cluster Energy [MeV]");
197
198 fhClusterEnergyMomentVSEnergy =
199 new TH2D("ClusterEnergyMomentVSEnergy", "Cluster Energy Moment VS Energy", 100, 0, 100, 60, 0, 600);
200 fhClusterEnergyMomentVSEnergy->GetXaxis()->SetTitle("|#vec{EM}| [cm]");
201 fhClusterEnergyMomentVSEnergy->GetYaxis()->SetTitle("Cluster Energy [MeV]");
202
203 fhClusterEnergyMoment = new TH1D("ClusterEnergyMoment", "Cluster Energy Moment", 100, 0, 100);
204 fhClusterMaxEnergyDigiMinusFirstDigiMag =
205 new TH1D("ClusterMaxEnergyDigiMinusFirstDigiMag", "ClusterMaxEnergyDigiMinusFirstDigiMag", 100, 0, 1000);
206
207 fhClusterEnergyMomentVSClusterSize =
208 new TH2D("ClusterEnergyMomentVSClusterSize", "Cluster Energy Moment VS Cluster Size", 100, 0, 100, 30, 0, 30);
209
210 fhEToFVSEelastic = new TH2D("EToFVSEelastic", "EToFVSEelastic", 50, 0, 1000, 50, 10, 1000);
211 fhScatteredNEnergyVSAngle = new TH2D("ScatteredNEnergyVSAngle", "ScatteredNEnergyVSAngle", 50, 0, 1000, 50, 0, 3.2);
212 fhScatteredNEnergyVSEdep = new TH2D("ScatteredNEnergyVSEdep", "ScatteredNEnergyVSEdep", 100, 0, 1000, 100, 0, 1000);
213
214 fhScatterAngleVSRecoilAngle =
215 new TH2D("ScatterAngleVSRecoilAngle", "ScatterAngleVSRecoilAngle", 50, 0, 3.2, 50, 0, 3.2);
216
217 fhSumAngleVSRatioErecoEtof = new TH2D("SumAngleVSRatioErecoEtof", "SumAngleVSRatioErecoEtof", 50, 0, 3.2, 50, 0, 2);
218 fhClusterEnergyVSScatteredRecoilAngle =
219 new TH2D("ClusterEnergyVSScatteredRecoilAngle", "ClusterEnergyVSScatteredRecoilAngle", 50, 0, 1000, 50, 0, 3.2);
220 fhClusterEnergyVSScatteredNeutronAngle = new TH2D(
221 "ClusterEnergyVSScatteredNeutronAngle", "ClusterEnergyVSScatteredNeutronAngle", 50, 0, 1000, 50, 0, 3.2);
222
223 fhElasticTargetMass = new TH2D("ElasticTargetMass", "ElasticTargetMass", 200, -50000, 50000, 100, 0, 1000);
224
225 /*fhtheh = new TH2D("theh", "", 100, 0, 1000, 100, 0, 2000);
226 ediff = new TH2D("ediff", "Reconstructed Proton Energy VS Measured Cluster Energy", 100, 0, 300, 100, 0, 300);
227 ediff->GetXaxis()->SetTitle("E_{p'}");
228 ediff->GetYaxis()->SetTitle("E_{cluster}");
229
230 etheta = new TH2D("etheta", "", 100, 0, 1000, 100, -1, 1);
231
232 fhProtonTrack*/
233
234 fhZ = new TH1D("Z", "Z", 5000, 0., 5000.);
235 fhZVSEToF = new TH2D("ZVSEToF", "Z vs EToF", 4000, 1000, 5000, 2000, 0, 2000);
236 fhDistFromCenter = new TH1D("DistFromCenter", "DistFromCenter", 500, 0., 500.);
237 fhDistFromCenterVSEToF = new TH2D("DistFromCenterVSEToF", "DistFromCenter vs EToF", 500, 0, 500, 2000, 0, 2000);
238
239 fhDeltaT = new TH1D("fhDeltaT", "T_{Last Digi} - T_{First Digi}", 1000, 0, 100);
240 fhForemostMinusFirstDigiTime =
241 new TH1D("fhForemostMinusFirstDigiTime", "T_{Foremost} - T_{First Digi}", 1000, -10, 10);
242
243 fhThetaEDigi = new TH2D("fhThetaEDigi", "fhThetaEDigi", 1000, -500, 500, 400, 0, 400);
244 fhThetaEDigiCosTheta = new TH2D("fhThetaEDigiCosTheta", "fhThetaEDigiCosTheta", 1000, -500, 500, 400, 0, 400);
245
246 hT = new TH1D("hT", "Cluster Digi Delta T", 3000, 0., 15.);
247 hTNeigh = new TH1D("hTNeigh", "Cluster Digi Neigh Delta T", 3000, 0., 15.);
248
249 return kSUCCESS;
250}
251
253{
254 fNeulandClustersBuffer = fNeulandClusters.get();
255 fNeulandClustersBuffer.erase(std::remove_if(fNeulandClustersBuffer.begin(),
256 fNeulandClustersBuffer.end(),
257 [&](R3BNeulandCluster& c) { return !(fClusterFilters.IsValid(&c)); }),
258 fNeulandClustersBuffer.end());
259
260 const auto nClusters = fNeulandClustersBuffer.size();
261
262 if (fIs3DTrackEnabled)
263 {
264 fh3->Reset("ICES");
265 for (const auto& cluster : fNeulandClustersBuffer)
266 {
267 const auto start = cluster.GetFirstHit().GetPosition();
268 // XYZ -> ZXY (side view)
269 fh3->Fill(start.Z(), start.X(), start.Y(), cluster.GetE());
270 }
271 }
272
273 const Double_t etot = std::accumulate(fNeulandClustersBuffer.begin(),
274 fNeulandClustersBuffer.end(),
275 0.,
276 [](Double_t sum, const R3BNeulandCluster& c) { return sum + c.GetE(); });
277
278 fhClusterNumberVSEnergy->Fill(etot, nClusters);
279
280 fhClusters->Fill(nClusters);
281 for (const auto& cluster : fNeulandClustersBuffer)
282 {
283 fhClusterTime->Fill(cluster.GetT());
284 fhClusterSize->Fill(cluster.GetSize());
285 fhClusterEnergy->Fill(cluster.GetE());
286 fhClusterRValue->Fill(std::log10(cluster.GetRCluster(fBeta)));
287 fhClusterEnergyVSSize->Fill(cluster.GetE(), cluster.GetSize());
288 fhClusterEToF->Fill(cluster.GetFirstHit().GetEToF());
289 fhClusterEToFVSEnergy->Fill(cluster.GetFirstHit().GetEToF(), cluster.GetE());
290 fhClusterEToFVSTime->Fill(cluster.GetFirstHit().GetEToF(), cluster.GetT());
291 fhClusterEVSTime->Fill(cluster.GetFirstHit().GetE(), cluster.GetT());
292
293 fhClusterEnergyVSEToF->Fill(cluster.GetE(), cluster.GetFirstHit().GetEToF());
294 fhClusterSizeVSEToF->Fill(cluster.GetSize(), cluster.GetFirstHit().GetEToF());
295 fhClusterEnergyVSSizeVSEToF->Fill(cluster.GetE(), cluster.GetSize(), cluster.GetFirstHit().GetEToF());
296 if (cluster.GetSize() > 2)
297 {
298 fhClusterForemostMinusCentroidVSEnergy->Fill(
299 (cluster.GetForemostHit().GetPosition() - cluster.GetEnergyCentroid()).Mag(), cluster.GetE());
300
301 fhClusterForemostMinusMaxEnergyDigiPosVSEnergy->Fill(
302 (cluster.GetForemostHit().GetPosition() - cluster.GetMaxEnergyHit().GetPosition()).Mag(),
303 cluster.GetE());
304
305 fhClusterCentroidMinusFirstDigiPosVSEnergy->Fill(
306 (cluster.GetEnergyCentroid() - cluster.GetFirstHit().GetPosition()).Mag(), cluster.GetE());
307
308 fhClusterMaxEnergyDigiMinusFirstDigiPosVSEnergy->Fill(
309 (cluster.GetMaxEnergyHit().GetPosition() - cluster.GetFirstHit().GetPosition()).Mag(), cluster.GetE());
310 fhClusterMaxEnergyDigiMinusCentroidVSEnergy->Fill(
311 (cluster.GetMaxEnergyHit().GetPosition() - cluster.GetEnergyCentroid()).Mag(), cluster.GetE());
312 fhClusterEnergyMomentVSEnergy->Fill(cluster.GetEnergyMoment(), cluster.GetE());
313 fhClusterEnergyMomentVSClusterSize->Fill(cluster.GetEnergyMoment(), cluster.GetSize());
314
315 fhClusterLastMinusFirstDigiMagVSEnergy->Fill(
316 (cluster.GetLastHit().GetPosition() - cluster.GetFirstHit().GetPosition()).Mag(), cluster.GetE());
317
318 fhClusterEnergyMoment->Fill(cluster.GetEnergyMoment());
319 fhClusterMaxEnergyDigiMinusFirstDigiMag->Fill(
320 (cluster.GetMaxEnergyHit().GetPosition() - cluster.GetFirstHit().GetPosition()).Mag());
321 }
322
323 fhZ->Fill(cluster.GetFirstHit().GetPosition().Z());
324 fhZVSEToF->Fill(cluster.GetFirstHit().GetPosition().Z(), cluster.GetEToF());
325 fhDistFromCenter->Fill(std::sqrt(std::pow(cluster.GetFirstHit().GetPosition().X(), 2) +
326 std::pow(cluster.GetFirstHit().GetPosition().Y(), 2)));
327 fhDistFromCenterVSEToF->Fill(std::sqrt(std::pow(cluster.GetFirstHit().GetPosition().X(), 2) +
328 std::pow(cluster.GetFirstHit().GetPosition().Y(), 2)),
329 cluster.GetEToF());
330 fhDeltaT->Fill(cluster.GetLastHit().GetT() - cluster.GetFirstHit().GetT());
331
332 fhForemostMinusFirstDigiTime->Fill(cluster.GetForemostHit().GetT() - cluster.GetFirstHit().GetT());
333
334 if (cluster.GetSize() > 4)
335 {
336 const auto theta = GetTheta(cluster);
337 for (const auto& digi : cluster.GetHits())
338 {
339 fhThetaEDigi->Fill(theta, digi.GetE());
340 fhThetaEDigiCosTheta->Fill(theta, digi.GetE() * std::cos(theta / rad2deg));
341 }
342 }
343 }
344
345 std::sort(fNeulandClustersBuffer.begin(),
346 fNeulandClustersBuffer.end(),
347 [](const R3BNeulandCluster& left, const R3BNeulandCluster& right) { return left.GetT() < right.GetT(); });
348
349 for (const auto& cluster : fNeulandClustersBuffer)
350 {
351 if (cluster.GetSize() >= 3)
352 {
353 fhENFromScatterVSEToF->Fill(Neuland::NeutronEnergyFromElasticProtonScattering(cluster), cluster.GetEToF());
354
355 fhClusterEnergyVSScatteredRecoilAngle->Fill(cluster.GetE(),
356 std::acos(Neuland::RecoilScatteringAngle(cluster)));
357 }
358 }
359
360 for (auto ita = fNeulandClustersBuffer.cbegin(); ita != fNeulandClustersBuffer.cend(); ita++)
361 {
362 for (auto itb = ita + 1; itb != fNeulandClustersBuffer.cend(); itb++)
363 {
364
365 fhScatteredNEnergyVSAngle->Fill(Neuland::ScatteredNeutronEnergy(*ita, *itb),
366 std::acos(Neuland::ScatteredNeutronAngle(*ita, *itb)));
367
368 fhScatteredNEnergyVSEdep->Fill(Neuland::ScatteredNeutronEnergy(*ita, *itb), (*ita).GetE());
369
370 const Double_t EelasticHeavy = Neuland::NeutronEnergyFromElasticScattering(*ita, *itb, 11000);
371 fhEToFVSEelastic->Fill((*ita).GetFirstHit().GetEToF(), EelasticHeavy);
372 // const Double_t EelasticProton = Neuland::NeutronEnergyFromElasticScattering(*ita, *itb, 1000);
373 // fhEToFVSEelastic->Fill((*ita)->GetFirstHit().GetEToF(), EelasticProton);
374
375 if (Neuland::ScatteredNeutronEnergy(*ita, *itb) > 10.)
376 {
377 fhElasticTargetMass->Fill(Neuland::ElasticScatteringTargetMass(*ita, *itb), (*ita).GetE());
378 fhClusterEnergyVSScatteredNeutronAngle->Fill((*ita).GetE(),
379 std::acos(Neuland::ScatteredNeutronAngle(*ita, *itb)));
380
381 if ((*ita).GetSize() >= 3)
382 {
383 fhScatterAngleVSRecoilAngle->Fill(std::acos(Neuland::ScatteredNeutronAngle(*ita, *itb)),
384 std::acos(Neuland::RecoilScatteringAngle(*ita)));
385 fhSumAngleVSRatioErecoEtof->Fill(
386 std::acos(Neuland::ScatteredNeutronAngle(*ita, *itb)) +
387 std::acos(Neuland::RecoilScatteringAngle(*ita)),
388 Neuland::NeutronEnergyFromElasticProtonScattering(*ita) / (*ita).GetFirstHit().GetEToF());
389 }
390 }
391 }
392 }
393
394 for (const auto& cluster : fNeulandClustersBuffer)
395 {
396 const auto& digis = cluster.GetHits();
397 for (auto it1 = digis.begin(); it1 != digis.end(); it1++)
398 {
399 for (auto it2 = it1 + 1; it2 != digis.end(); it2++)
400 {
401 if (std::abs(it1->GetPosition().X() - it2->GetPosition().X()) < 7.5 &&
402 std::abs(it1->GetPosition().Y() - it2->GetPosition().Y()) < 7.5 &&
403 std::abs(it1->GetPosition().Z() - it2->GetPosition().Z()) < 7.5)
404 {
405 hTNeigh->Fill(std::abs(it1->GetT() - it2->GetT()));
406 }
407 hT->Fill(std::abs(it1->GetT() - it2->GetT()));
408 }
409 }
410 }
411}
412
414{
415 TDirectory* tmp = gDirectory;
416 FairRootManager::Instance()->GetOutFile()->cd();
417
418 gDirectory->mkdir(fOutput);
419 gDirectory->cd(fOutput);
420
421 fhClusters->Write();
422 fhClusterSize->Write();
423 fhClusterEnergyVSSize->Write();
424 fhClusterEnergy->Write();
425 fhClusterRValue->Write();
426 fhClusterTime->Write();
427 fhClusterNumberVSEnergy->Write();
428 fhClusterEToF->Write();
429 fhClusterEToFVSEnergy->Write();
430 fhClusterEToFVSTime->Write();
431 fhElasticTargetMass->Write();
432
433 fhClusterSizeVSEToF->Write();
434 fhClusterEnergyVSEToF->Write();
435 // fhClusterEnergyVSSizeVSEToF->Write();
436
437 fhClusterForemostMinusCentroidVSEnergy->Write();
438 fhClusterForemostMinusMaxEnergyDigiPosVSEnergy->Write();
439 fhClusterCentroidMinusFirstDigiPosVSEnergy->Write();
440 fhClusterMaxEnergyDigiMinusFirstDigiPosVSEnergy->Write();
441 fhClusterMaxEnergyDigiMinusCentroidVSEnergy->Write();
442 fhClusterEnergyMomentVSEnergy->Write();
443 fhClusterEnergyMomentVSClusterSize->Write();
444 fhClusterEVSTime->Write();
445 fhClusterEnergyMoment->Write();
446 fhClusterMaxEnergyDigiMinusFirstDigiMag->Write();
447 fhClusterLastMinusFirstDigiMagVSEnergy->Write();
448
449 fhENFromScatterVSEToF->Write();
450
451 fhEToFVSEelastic->Write();
452 fhScatteredNEnergyVSAngle->Write();
453 fhScatteredNEnergyVSEdep->Write();
454 fhScatterAngleVSRecoilAngle->Write();
455 fhSumAngleVSRatioErecoEtof->Write();
456 fhClusterEnergyVSScatteredRecoilAngle->Write();
457 fhClusterEnergyVSScatteredNeutronAngle->Write();
458
459 fhZ->Write();
460 fhZVSEToF->Write();
461 fhDistFromCenter->Write();
462 fhDistFromCenterVSEToF->Write();
463 fhDeltaT->Write();
464 fhForemostMinusFirstDigiTime->Write();
465
466 fhThetaEDigi->Write();
467 fhThetaEDigiCosTheta->Write();
468
469 hT->Write();
470 hTNeigh->Write();
471
472 gDirectory = tmp;
473}
474
ClassImp(R3B::Neuland::Cal2HitPar)
Double_t GetTheta(const R3BMCTrack *mcTrack)
static const Double_t c
void Exec(Option_t *) override
R3BNeulandClusterMon(TString input="NeulandClusters", TString output="NeulandClusterMon", const Option_t *option="")
auto Init() -> InitStatus override
auto NeutronEnergyFromElasticProtonScattering(const R3BNeulandCluster &cluster) -> double
auto ElasticScatteringTargetMass(const R3BNeulandCluster &first, const R3BNeulandCluster &second) -> double
auto NeutronEnergyFromElasticScattering(const R3BNeulandCluster &first, const R3BNeulandCluster &second, double target_mass) -> double
Calculate neutron kinetic energy from the elastic scattering.
auto ScatteredNeutronEnergy(const R3BNeulandCluster &first, const R3BNeulandCluster &second) -> double
Calculate neutron kinetic energy from the scattering.
auto ScatteredNeutronAngle(const R3BNeulandCluster &first, const R3BNeulandCluster &second) -> double
Calculate the scattering angle of the neutron.
auto RecoilScatteringAngle(const R3BNeulandCluster &cluster) -> double