Skip to content

Commit 7fd1ea1

Browse files
Move FairModule::fAllSensitiveVolumes to FairRunSim
Instead of using a thread local, let's use something that is tied to the "simulaion session" (as noted in some places). FairRunSim has all the Modules anyway, so let it also have the list of sensitive Volumes.
1 parent 54ad232 commit 7fd1ea1

File tree

5 files changed

+18
-13
lines changed

5 files changed

+18
-13
lines changed

fairroot/base/sim/FairDetector.cxx

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -14,6 +14,7 @@
1414

1515
#include "FairGeoNode.h" // for FairGeoNode
1616
#include "FairRootManager.h"
17+
#include "FairRunSim.h"
1718
#include "FairVolume.h" // for FairVolume
1819

1920
#include <TFolder.h> // for TFolder
@@ -97,7 +98,7 @@ void FairDetector::Initialize()
9798
FairGeoNode* fN;
9899
TString cutName;
99100
TString copysign = "#";
100-
for (auto aVol : fAllSensitiveVolumes) {
101+
for (auto aVol : FairRunSim::Instance()->fAllSensitiveVolumes) {
101102
cutName = aVol->GetName();
102103
Ssiz_t pos = cutName.Index(copysign, 1);
103104
if (pos > 1) {

fairroot/base/sim/FairMCApplication.cxx

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -886,7 +886,7 @@ void FairMCApplication::InitGeometry()
886886
fMCEventHeader->SetRunID(runId);
887887

888888
// Fill sensitive volumes in fVolMap
889-
for (auto fv : FairModule::fAllSensitiveVolumes) {
889+
for (auto fv : fRun->fAllSensitiveVolumes) {
890890
if (!fv) {
891891
continue;
892892
}

fairroot/base/sim/FairModule.cxx

Lines changed: 5 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -20,7 +20,7 @@
2020
#include "FairGeoNode.h" // for FairGeoNode
2121
#include "FairGeoParSet.h" // for FairBaseParSet
2222
#include "FairMCApplication.h"
23-
#include "FairRun.h" // for FairRun
23+
#include "FairRunSim.h" // for FairRun, FairRunSim
2424
#include "FairRuntimeDb.h" // for FairRuntimeDb
2525
#include "FairVolume.h" // for FairVolume
2626
#include "FairVolumeList.h" // for FairVolumeList
@@ -47,8 +47,6 @@
4747
#include <map>
4848
#include <memory>
4949

50-
thread_local std::vector<FairVolume*> FairModule::fAllSensitiveVolumes;
51-
5250
void FairModule::ConstructGeometry()
5351
{
5452
LOG(warn)
@@ -214,11 +212,12 @@ void FairModule::ProcessNodes(TList* aList)
214212
vList = new FairVolumeList();
215213
}
216214

215+
FairRunSim& run = *(FairRunSim::Instance());
217216
TListIter iter(aList);
218217
FairGeoNode* node = nullptr;
219218
FairGeoNode* MotherNode = nullptr;
220219
FairVolume* volume = nullptr;
221-
FairRuntimeDb* rtdb = FairRun::Instance()->GetRuntimeDb();
220+
FairRuntimeDb* rtdb = run.GetRuntimeDb();
222221
FairGeoParSet* par = static_cast<FairGeoParSet*>(rtdb->getContainer("FairGeoParSet"));
223222
TObjArray* fNodes = par->GetGeoNodes();
224223
while ((node = static_cast<FairGeoNode*>(iter.Next()))) {
@@ -240,7 +239,7 @@ void FairModule::ProcessNodes(TList* aList)
240239
if (node->isSensitive() && fActive) {
241240
volume->setModId(fModId);
242241
volume->SetModule(this);
243-
fAllSensitiveVolumes.push_back(volume);
242+
run.fAllSensitiveVolumes.push_back(volume);
244243
aVol = dynamic_cast<FairGeoVolume*>(node);
245244
fNodes->AddLast(aVol);
246245
fNbOfSensitiveVol++;
@@ -260,7 +259,7 @@ void FairModule::AddSensitiveVolume(TGeoVolume* v)
260259
vList->addVolume(volume);
261260
volume->setModId(fModId);
262261
volume->SetModule(this);
263-
fAllSensitiveVolumes.push_back(volume);
262+
FairRunSim::Instance()->fAllSensitiveVolumes.push_back(volume);
264263
fNbOfSensitiveVol++;
265264
}
266265
}

fairroot/base/sim/FairModule.h

Lines changed: 0 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -22,11 +22,9 @@
2222
#include <TList.h> // for TList (ptr only), TListIter
2323
#include <TNamed.h> // for TNamed
2424
#include <TObjArray.h> // for TObjArray
25-
#include <TRefArray.h> // for TRefArray
2625
#include <TString.h> // for TString, operator!=
2726
#include <TVirtualMC.h>
2827
#include <string>
29-
#include <vector>
3028

3129
class FairVolumeList;
3230
class FairVolume;
@@ -139,8 +137,6 @@ class FairModule : public TNamed
139137
static thread_local inline FairVolumeList* vList{nullptr}; //!
140138
/**total number of volumes in a simulaion session*/
141139
static thread_local inline Int_t fNbOfVolumes{0}; //!
142-
/**list of all sensitive volumes in a simulaion session*/
143-
static thread_local std::vector<FairVolume*> fAllSensitiveVolumes; //!
144140

145141
TString fMotherVolumeName{""}; //!
146142
FairVolume* getFairVolume(FairGeoNode* fNode);

fairroot/base/steer/FairRunSim.h

Lines changed: 10 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,5 +1,5 @@
11
/********************************************************************************
2-
* Copyright (C) 2014-2023 GSI Helmholtzzentrum fuer Schwerionenforschung GmbH *
2+
* Copyright (C) 2014-2024 GSI Helmholtzzentrum fuer Schwerionenforschung GmbH *
33
* *
44
* This software is distributed under the terms of the *
55
* GNU Lesser General Public Licence (LGPL) version 3, *
@@ -21,12 +21,14 @@
2121
#include <functional>
2222
#include <memory>
2323
#include <utility>
24+
#include <vector>
2425

2526
class FairField;
2627
class FairMCEventHeader;
2728
class FairMesh;
2829
class FairModule;
2930
class FairPrimaryGenerator;
31+
class FairVolume;
3032

3133
/**
3234
* \brief Configure the Simulation session
@@ -37,6 +39,10 @@ class FairPrimaryGenerator;
3739
*/
3840
class FairRunSim : public FairRun
3941
{
42+
friend class FairModule;
43+
friend class FairDetector;
44+
friend class FairMCApplication;
45+
4046
public:
4147
/** default ctor*/
4248
FairRunSim(Bool_t isMaster = kTRUE);
@@ -223,6 +229,9 @@ class FairRunSim : public FairRun
223229

224230
bool fWasMT{false}; //! /** Actual MT mode used */
225231

232+
/**list of all sensitive volumes in a simulaion session*/
233+
std::vector<FairVolume*> fAllSensitiveVolumes; //!
234+
226235
protected:
227236
Int_t count{0}; //!< Internal counter
228237
FairMCApplication* fApp{nullptr}; //!< Main VMC application

0 commit comments

Comments
 (0)