KaliVeda
Toolkit for HIC analysis
KVFAZIARawDataReconstructor.cpp
1 //Created by KVClassFactory on Fri Jan 23 14:17:23 2015
2 //Author: ,,,
3 
4 #include "KVFAZIARawDataReconstructor.h"
5 #include "KVDataAnalyser.h"
6 #include "KVDataSet.h"
7 #include "KVFAZIADB.h"
8 #include "KVDataRepositoryManager.h"
9 #include "KVDataRepository.h"
10 #include "RVersion.h"
11 #include "KVFAZIADetector.h"
12 #include "KVFAZIA.h"
13 #include "KVSignal.h"
14 #include "TSystem.h"
15 
17 
18 // BEGIN_HTML <!--
20 /* -->
21 <h2>KVFAZIARawDataReconstructor</h2>
22 <h4>Handle reconstruction of FAZIA events</h4>
23 <!-- */
24 // --> END_HTML
26 
27 
28 
32  : taskname("Reconstruction"), datatype("recon")
33 {
34  //Default constructor
35  file = 0;
36  tree = 0;
37  recev = 0;
38  nb_recon = 0;
39 }
40 
41 
42 
45 
47 {
48  // Destructor
50 }
51 
52 
53 
54 
68 
70 {
71  // Creates new ROOT file with TTree for reconstructed events.
72  // By default this file will be written in the same data repository as the raw data file we are reading.
73  // This can be changed by setting the environment variable(s):
74  //
75  // Reconstruction.DataAnalysisTask.OutputRepository: [name of repository]
76  // [name of dataset].Reconstruction.DataAnalysisTask.OutputRepository: [name of repository]
77  //
78  // If no value is set for the current dataset (second variable), the value of the
79  // first variable will be used. If neither is defined, the new file will be written in the same repository as
80  // the raw file (if possible, i.e. if repository is not remote).
81 
82  // Create new KVReconstructedEvent filled with KVFAZIAReconNuc object
83  // used to reconstruct & store events
84 
85  if (!recev) recev = new KVReconstructedEvent(50);
86 
87  // get dataset to which we must associate new run
88  KVDataSet* OutputDataset =
89  gDataRepositoryManager->GetDataSet(gDataSet->GetOutputRepository(taskname), gDataSet->GetName());
90 
91  file = OutputDataset->NewRunfile(datatype.Data(), GetCurrentRunNumber());
92 
93  std::cout << "Writing \"" << datatype.Data() << "\" events in ROOT file " << file->GetName() << std::endl;
94 
95  //tree for reconstructed events
96  tree = new TTree("ReconstructedEvents", Form("%s : %s : %s events created from raw data",
97  gFaziaDB->GetRun(GetCurrentRunNumber())->GetName(),
98  gFaziaDB->GetRun(GetCurrentRunNumber())->GetTitle(),
99  datatype.Data())
100  );
101 
102  //leaves for reconstructed events
103  KVEvent::MakeEventBranch(tree, "FAZIAReconEvent", recev);
104 
105  Info("InitRun", "Created reconstructed data tree %s : %s", tree->GetName(), tree->GetTitle());
106  nb_recon = 0;
107 
108 }
109 
110 
111 
112 
120 
122 {
123  // Analysis of event.
124  // RawData TTree is filled with values of parameters for event.
125  // If event has INDRA trigger information (IsIndraEvent()==kTRUE)
126  // then
127  // *) event reconstruction is performed for 'Physics' events
128  // *) or the GeneTree is filled with pulser/laser data for 'Gene' events
129 
130 
131  recev->SetNumber(GetEvent()->GetNumber());
132  gFazia->ReconstructEvent(recev, GetDetectorEvent());
133 
134  ExtraProcessing();
135 
136  nb_recon++;
137  tree->Fill();
138 
139  recev->Clear();
140 
141  return kTRUE;
142 }
143 
144 
145 
146 
148 
150 {
151  KVString label = "";
152 
153  KVFAZIADetector* det = 0;
154  KVSignal* sig = 0;
155  KVReconstructedNucleus* recnuc = 0;
156  while ((recnuc = recev->GetNextParticle())) {
158  KVGeoDetectorNode* _node;
159  while ((_node = recnuc->GetReconstructionTrajectory()->GetNextNode())) {
160  auto det = dynamic_cast<KVFAZIADetector*>(_node->GetDetector());
161  TIter next_s(det->GetListOfSignals());
162  while ((sig = (KVSignal*)next_s())) {
163  if (sig->HasFPGA()) {
164  for (Int_t ii = 0; ii < sig->GetNFPGAValues(); ii += 1) {
165  //SI2-T3-Q1-B003.Q2.RawAmplitude=14
166  if (ii == 0) label = "FPGAEnergy";
167  if (ii == 1) label = "FastFPGAEnergy"; //only for CsI Q3
168  TString ene = GetEvent()->GetFPGAEnergy(
169  det->GetBlockNumber(),
170  det->GetQuartetNumber(),
171  det->GetTelescopeNumber(),
172  sig->GetType(),
173  ii
174  );
175 
176  recnuc->GetParameters()->SetValue(
177  Form("%s.%s.%s", det->GetName(), sig->GetName(), label.Data()),
178  ene.Data()
179  );
180  }
181  }
182  if (!sig->PSAHasBeenComputed()) {
183  sig->TreateSignal();
184  }
185 
186  sig->GetPSAResult(det);
187 // if (psa.GetNpar()) {
188 // *(recnuc->GetParameters()) += psa;
189 // }
190  }
191  }
192  }
193 }
194 
195 
196 
197 
199 
201 {
202  SafeDelete(recev);
203 
204  std::cout << std::endl << " *** Number of reconstructed FAZIA events : "
205  << nb_recon << " ***" << std::endl << std::endl;
206  file->cd();
207  gDataAnalyser->WriteBatchInfo(tree);
208  tree->Write();//write tree to file
209 
210  // get dataset to which we must associate new run
211  KVDataSet* OutputDataset =
212  gDataRepositoryManager->GetDataSet(gDataSet->GetOutputRepository(taskname), gDataSet->GetName());
213  //add new file to repository
214  OutputDataset->CommitRunfile(datatype.Data(), GetCurrentRunNumber(), file);
215 
216  ProcInfo_t pid;
217  if (gSystem->GetProcInfo(&pid) == 0) {
218  std::cout << " ------------- Process infos -------------" << std::endl;
219  printf(" CpuSys = %f s. CpuUser = %f s. ResMem = %f MB VirtMem = %f MB\n",
220  pid.fCpuSys, pid.fCpuUser, pid.fMemResident / 1024., pid.fMemVirtual / 1024.);
221  }
222 
223 }
224 
225 
int Int_t
#define SafeDelete(p)
bool Bool_t
constexpr Bool_t kTRUE
char * Form(const char *fmt,...)
R__EXTERN TSystem * gSystem
virtual void SetNumber(UInt_t num)
Definition: KVBase.h:215
void WriteBatchInfo(TTree *)
KVDataSet * GetDataSet(const Char_t *repository, const Char_t *dataset) const
Return pointer to named dataset in the given repository.
Manage an experimental dataset corresponding to a given experiment or campaign.
Definition: KVDataSet.h:36
TString GetOutputRepository(const Char_t *taskname) const
Definition: KVDataSet.cpp:2086
void CommitRunfile(const KVString &type, const run_index_t &run, TFile *file)
Definition: KVDataSet.cpp:1372
TFile * NewRunfile(const KVString &type, const run_index_t &run)
Definition: KVDataSet.cpp:1112
static void MakeEventBranch(TTree *tree, const TString &branchname, T &event, Int_t bufsize=10000000)
Definition: KVEvent.h:210
void Clear(Option_t *opt="") override
Definition: KVEvent.h:238
KVFAZIADBRun * GetRun(Int_t run) const
Definition: KVFAZIADB.h:48
Base class for FAZIA detectors.
const KVSeqCollection * GetListOfSignals() const
Int_t GetQuartetNumber() const
Int_t GetBlockNumber() const
Int_t GetTelescopeNumber() const
Reconstruction of raw FAZIA data.
virtual ~KVFAZIARawDataReconstructor()
Destructor.
Int_t nb_recon
number of reconstructed events
const Char_t * GetFPGAEnergy(Int_t blk, Int_t qua, Int_t tel, TString signaltype, Int_t idx=0)
Int_t GetCurrentRunNumber() const
Definition: KVFAZIAReader.h:74
KVFAZIARawEvent * GetEvent() const
Definition: KVFAZIAReader.h:82
KVDetectorEvent * GetDetectorEvent() const
Definition: KVFAZIAReader.h:86
KVGeoDetectorNode * GetNextNode() const
void IterateFrom(const KVGeoDetectorNode *node0=nullptr) const
Information on relative positions of detectors & particle trajectories.
KVDetector * GetDetector() const
void SetValue(const Char_t *name, value_type value)
KVNameValueList * GetParameters() const
Definition: KVParticle.h:818
Event containing KVReconstructedNucleus nuclei reconstructed from hits in detectors.
Nuclei reconstructed from data measured by a detector array .
const KVReconNucTrajectory * GetReconstructionTrajectory() const
const Char_t * GetType() const
Definition: KVSignal.h:98
virtual void TreateSignal()
Definition: KVSignal.cpp:406
virtual void GetPSAResult(KVDetector *) const
Definition: KVSignal.h:155
Int_t GetNFPGAValues() const
Definition: KVSignal.h:129
Bool_t PSAHasBeenComputed() const
Definition: KVSignal.h:159
Bool_t HasFPGA() const
Definition: KVSignal.h:134
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:73
Particle * GetNextParticle(Option_t *opt="") const
const char * GetName() const override
virtual const char * GetName() const
virtual const char * GetTitle() const
virtual void Info(const char *method, const char *msgfmt,...) const
const char * Data() const
virtual int GetProcInfo(ProcInfo_t *info) const
Long_t fMemVirtual
Float_t fCpuSys
Long_t fMemResident
Float_t fCpuUser
ClassImp(TPyArg)