R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeuland.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
14#include "R3BNeuland.h"
15#include "FairRun.h"
16#include "FairRuntimeDb.h"
17#include "R3BDetectorList.h"
18#include "R3BException.h"
19#include "R3BMCStack.h"
20#include "R3BNeulandGeoPar.h"
21#include "R3BNeulandPoint.h"
22#include "TClonesArray.h"
23#include "TGeoManager.h"
24#include "TVirtualMC.h"
25#include <FairGeoLoader.h>
26#include <FairRootManager.h>
27#include <FairVolume.h>
28#include <R3BLogger.h>
29#include <Rtypes.h>
30#include <RtypesCore.h>
31#include <TCollection.h>
32#include <TGeoMatrix.h>
33#include <TGeoNode.h>
34#include <TMathBase.h>
35#include <TString.h>
36#include <fairlogger/Logger.h>
37#include <fmt/core.h>
38#include <memory>
39#include <string>
40
41// Initialize variables from Birk' s Law
42constexpr auto seconds_to_nanoseconds = 1e9;
43constexpr auto BirkdP = 1.032;
44constexpr auto BirkC1 = 0.013 / BirkdP;
45constexpr auto BirkC2 = 9.6e-6 / (BirkdP * BirkdP);
46
47namespace
48{
49 inline auto GetLightYield(const int charge, const double length, const double edep) -> double
50 {
51 // Apply Birk's law ( Adapted from G3BIRK/Geant3)
52 if (charge != 0 && length > 0)
53 {
54 auto birkC1Mod = BirkC1;
55
56 // Apply correction for higher charge states
57 if (TMath::Abs(charge) >= 2)
58 {
59 birkC1Mod *= 7.2 / 12.6; // NOLINT
60 }
61
62 const double dedxcm = 1000. * edep / length;
63 const double lightYield = edep / (1. + birkC1Mod * dedxcm + BirkC2 * dedxcm * dedxcm);
64 return lightYield;
65 }
66 return edep; // Rarely very small energy depositions have no length?
67 }
68} // namespace
69
74
75R3BNeuland::R3BNeuland(const TString& geoFile, const TGeoTranslation& trans, const TGeoRotation& rot)
76 : R3BNeuland(geoFile, { trans, rot })
77{
78}
79
80R3BNeuland::R3BNeuland(const TString& geoFile, const TGeoCombiTrans& combi)
81 : FairDetector{ "R3BNeuland", true, kNEULAND }
82 , rot_trans_{ combi }
83 , geo_file_{ geoFile.Data() }
84{
85}
86
87R3BNeuland::R3BNeuland(int nDP, const TGeoTranslation& trans, const TGeoRotation& rot)
88 : R3BNeuland(nDP, { trans, rot })
89{
90}
91
92R3BNeuland::R3BNeuland(const int nDP, const TGeoCombiTrans& combi)
93 : R3BNeuland(fmt::format("neuland_v3_{}dp.geo.root", nDP), combi)
94{
95 num_of_planes_ = 2 * nDP;
96}
97
99{
100 LOG(info) << "R3BNeuland initialization ...";
101
102 FairDetector::Initialize();
103
105 reset_values();
106}
107
108auto R3BNeuland::ProcessHits(FairVolume* /*v*/) -> bool
109{
110 // New hit in detector
111 if (gMC->IsTrackEntering())
112 {
114 {
115 LOG(warn) << "R3BNeuland: Incomplete hit discarded";
116 reset_values();
117 }
118
119 is_last_hit_done_ = kFALSE;
120 energy_loss_ = 0.;
121 light_yield_ = 0.;
122 time_ = gMC->TrackTime() * seconds_to_nanoseconds;
123 length_ = gMC->TrackLength();
124 gMC->TrackPosition(pos_in_);
125 gMC->TrackMomentum(mom_in_);
126 gMC->CurrentVolOffID(1, paddle_id_);
127
128 particle_id_ = gMC->TrackPid();
129 track_pid_map_.emplace(gMC->GetStack()->GetCurrentTrackNumber(), gMC->TrackPid());
130 if (auto search = track_pid_map_.find(gMC->GetStack()->GetCurrentParentTrackNumber());
131 search != track_pid_map_.end())
132 {
133 parent_particle_id_ = search->first;
134 }
135 }
136
137 // Sum energy loss for all steps in the active volume
138 energy_loss_ += gMC->Edep();
139 light_yield_ += GetLightYield(gMC->TrackCharge(), gMC->TrackStep(), gMC->Edep());
140
141 // Set additional parameters at exit of active volume. Create R3BNeulandPoint.
142 if (gMC->IsTrackExiting() || gMC->IsTrackStop() || gMC->IsTrackDisappeared())
143 {
144 // Do not save a hit if no energy deposited
145 constexpr auto minimum_energy_cutoff = 1e-20;
146 if (energy_loss_ < minimum_energy_cutoff || light_yield_ < minimum_energy_cutoff)
147 {
148 reset_values();
149 return kTRUE;
150 }
151
152 track_id_ = gMC->GetStack()->GetCurrentTrackNumber();
153 gMC->TrackPosition(pos_out_);
154 gMC->TrackMomentum(mom_out_);
155
156 // Add Point
157 LOGP(debug,
158 "R3BNeuland: Adding Point at (\"{}\", \"{}\", \"{}\") cm, paddle {}, track {}, energy loss {} GeV",
159 pos_in_.X(),
160 pos_in_.Y(),
161 pos_in_.Z(),
163 track_id_,
165 gMC->GetStack()->GetCurrentParentTrackNumber());
166 auto* neuland_point =
167 dynamic_cast<R3BNeulandPoint*>(tca_points_buffer_->ConstructedAt(tca_points_buffer_->GetEntriesFast()));
168 neuland_point->SetTrackID(track_id_);
169 neuland_point->SetDetectorID(paddle_id_);
170 neuland_point->SetPosition(pos_in_.Vect());
171 neuland_point->SetMomentum(mom_in_.Vect());
172 neuland_point->SetTime(time_);
173 neuland_point->SetLength(length_);
174 neuland_point->SetEnergyLoss(energy_loss_);
175 neuland_point->SetEventID(gMC->CurrentEvent());
176 neuland_point->SetLightYield(light_yield_);
177 neuland_point->SetParticleId(particle_id_);
178 neuland_point->SetParentParticleId(parent_particle_id_);
179
180 // Increment number of LandPoints for this track
181 auto* stack = dynamic_cast<R3BStack*>(gMC->GetStack());
182 stack->AddPoint(kNEULAND);
183 reset_values();
184 }
185
186 return kTRUE;
187}
188
189auto R3BNeuland::CheckIfSensitive(std::string name) -> bool { return name == "volBC408"; }
190
192{
193 if (fVerboseLevel != 0)
194 {
195 Print();
196 }
197 Reset();
198}
199
201{
202 auto& points = neuland_points_.get();
203 points.reserve(tca_points_buffer_->GetEntriesFast());
204 for (auto* point : TRangeDynCast<R3BNeulandPoint>(tca_points_buffer_.get()))
205 {
206 points.push_back(*point);
207 }
208}
209
210void R3BNeuland::Print(Option_t* /*unused*/) const
211{
212 LOG(info) << "R3BNeuland: " << neuland_points_.get_constref().size() << " Neuland Points registered in this event";
213}
214
216{
217 neuland_points_.clear();
218 tca_points_buffer_->Clear();
219 reset_values();
220 track_pid_map_.clear();
221}
222
224{
225 is_last_hit_done_ = kTRUE;
226 track_id_ = 0;
227 paddle_id_ = -1;
228 time_ = 0;
229 length_ = 0.;
230 energy_loss_ = 0.;
231 light_yield_ = 0.;
232 pos_in_.Clear();
233 pos_out_.Clear();
234 mom_in_.Clear();
235 mom_out_.Clear();
236}
237
239{
240 FairRuntimeDb* rtdb = FairRun::Instance()->GetRuntimeDb();
241 neuland_geo_par_ = dynamic_cast<R3BNeulandGeoPar*>(rtdb->getContainer("R3BNeulandGeoPar"));
242
243 // Really bad way to find the Neuland *node* (not the volume!)
244 TGeoNode* geoNodeNeuland = nullptr;
245 for (int i{}; i < gGeoManager->GetTopNode()->GetNdaughters(); i++)
246 {
247 if (TString(gGeoManager->GetTopNode()->GetDaughter(i)->GetVolume()->GetName()) == "volNeuland")
248 {
249 geoNodeNeuland = gGeoManager->GetTopNode()->GetDaughter(i);
250 break;
251 }
252 }
253
254 if (geoNodeNeuland == nullptr)
255 {
256 LOG(fatal) << "volNeuland not found";
257 }
258
259 neuland_geo_par_->SetNeulandGeoNode(geoNodeNeuland);
260 neuland_geo_par_->setChanged();
261}
262
263auto R3BNeuland::GetCollection(int iColl) const -> TClonesArray*
264{
265 if (iColl == 0)
266 {
267 return tca_points_buffer_.get();
268 }
269 return nullptr;
270}
271
273
275{
277 {
278 create_geo();
279 }
280 else
281 {
283 }
284}
285
287{
288
289 auto* geo_loader = FairGeoLoader::Instance();
290 if (geo_loader == nullptr)
291 {
292 geo_loader = std::make_unique<FairGeoLoader>("TGeo", "FairGeoLoader").release();
293 }
294 auto* neuland_geo = geo_creator_.construct_volume(num_of_planes_, geo_loader);
295 if (auto* top_volume = gGeoManager->GetTopVolume(); top_volume != nullptr)
296 {
297 // Use copy_id 0 since only one neuland is needed
298 auto* neuland_node = top_volume->AddNode(neuland_geo, 0, rot_trans_.MakeClone());
299
300 // This will set each bar as a sensitive volume
301 ExpandNode(neuland_node);
302 }
303 else
304 {
305 throw R3B::runtime_error("Top volume from gGeoManager is nullptr!");
306 }
307}
308
310{
311 SetGeometryFileName(geo_file_.c_str());
312 if (!GetGeometryFileName().EndsWith(".root"))
313 {
314 R3BLOG(fatal, GetName() << " (which is a " << ClassName() << ") geometry file is not specified");
315 }
316 R3BLOG(info,
317 fmt::format("Constructing {} (which is a {}) geometry from ROOT file {} ...",
318 GetName(),
319 ClassName(),
320 GetGeometryFileName().Data()));
321 ConstructRootGeometry();
322 if (not rot_trans_.IsIdentity())
323 {
324 if (auto* top_node = gGeoManager->GetTopNode(); top_node != nullptr)
325 {
326 auto* neuland_node = top_node->GetDaughter(gGeoManager->GetTopNode()->GetNdaughters() - 1);
327 dynamic_cast<TGeoNodeMatrix*>(neuland_node)->SetMatrix(rot_trans_.MakeClone());
328 }
329 else
330 {
331 throw R3B::runtime_error("Top node from gGeoManager is nullptr!");
332 }
333 }
334}
335
@ kNEULAND
#define R3BLOG(severity, x)
Definition R3BLogger.h:32
constexpr auto BirkdP
constexpr auto seconds_to_nanoseconds
constexpr auto BirkC2
ClassImp(R3BNeuland)
constexpr auto BirkC1
NeuLAND geometry parameter storage.
NeuLAND detector simulation class.
Definition R3BNeuland.h:49
bool is_last_hit_done_
Flag to check if last hit finished.
Definition R3BNeuland.h:98
auto CheckIfSensitive(std::string name) -> bool override
R3B::OutputVectorConnector< R3BNeulandPoint > neuland_points_
Output data written to the ROOT file.
Definition R3BNeuland.h:118
int track_id_
A buffer for the track IDs.
Definition R3BNeuland.h:102
TGeoCombiTrans rot_trans_
Rotation and tranlation of NeuLAND detector.
Definition R3BNeuland.h:124
TLorentzVector pos_out_
A buffer for the position where the track goes out of the volume.
Definition R3BNeuland.h:111
double energy_loss_
A buffer for the energy losses.
Definition R3BNeuland.h:108
auto ProcessHits(FairVolume *=nullptr) -> bool override
int paddle_id_
A buffer for the paddle ids.
Definition R3BNeuland.h:103
std::string geo_file_
Geometry file name. The file must be present if EnableAutoGeoBuild() is disabled.
Definition R3BNeuland.h:126
double time_
A buffer for time values.
Definition R3BNeuland.h:106
void reset_values()
void Register() override
TLorentzVector mom_out_
A buffer for the position when the track goes out of the volume.
Definition R3BNeuland.h:113
void FinishEvent() override
TLorentzVector pos_in_
A buffer for the position where the track goes into the volume.
Definition R3BNeuland.h:110
void EndOfEvent() override
bool is_geo_auto_built
Flag to check if geo-build is needed.
Definition R3BNeuland.h:97
R3B::Neuland::Geometry::Creator geo_creator_
Creating NeuLAND geometry if EnableAutoGeoBuild() is enabled.
Definition R3BNeuland.h:122
void create_geo()
void write_parameter_file()
void create_geo_from_root_file()
int particle_id_
A buffer for the particle IDs.
Definition R3BNeuland.h:104
R3BNeulandGeoPar * neuland_geo_par_
Output parameter containing the geometry information.
Definition R3BNeuland.h:120
int num_of_planes_
The number of planes.
Definition R3BNeuland.h:99
R3BNeuland()
Default constructor.
int parent_particle_id_
A buffer for the parent particle IDs.
Definition R3BNeuland.h:105
std::map< int, int > track_pid_map_
A map with the track ID as the key and the particle ID as the value.
Definition R3BNeuland.h:129
void ConstructGeometry() override
auto GetCollection(int iColl) const -> TClonesArray *override
void Reset() override
void Print(Option_t *="") const override
double length_
A buffer for track lengths.
Definition R3BNeuland.h:107
double light_yield_
A buffer for the light yields.
Definition R3BNeuland.h:109
void Initialize() override
TLorentzVector mom_in_
A buffer for the momentum when the track goes into the volume.
Definition R3BNeuland.h:112
std::unique_ptr< TClonesArray > tca_points_buffer_
The TCA data buffer used to be returned by virtual method GetCollection(int iColl).
Definition R3BNeuland.h:116
void SetLightYield(double light_yield)
R3BStack.h.
Definition R3BMCStack.h:55
void AddPoint(DetectorId iDet)
Increment number of points for the current track in a given detector.