FdEventLogger.cc
Go to the documentation of this file.
1 
14 #include <utl/Reader.h>
15 #include <utl/ErrorLogger.h>
16 #include <utl/AugerUnits.h>
17 #include <utl/MathConstants.h>
18 #include <utl/PhysicalConstants.h>
19 #include <utl/TabulatedFunction.h>
20 #include <utl/TabulatedFunctionErrors.h>
21 #include <utl/Trace.h>
22 #include <utl/MultiTabulatedFunction.h>
23 #include <utl/Point.h>
24 #include <utl/UTMPoint.h>
25 #include <utl/TimeInterval.h>
26 
27 #include <fwk/CentralConfig.h>
28 
29 #include <det/Detector.h>
30 #include <fdet/FDetector.h>
31 #include <fdet/Eye.h>
32 #include <fdet/Channel.h>
33 #include <fdet/Pixel.h>
34 
35 #include <evt/Event.h>
36 #include <evt/ShowerSimData.h>
37 #include <fevt/FEvent.h>
38 #include <fevt/Eye.h>
39 #include <fevt/EyeHeader.h>
40 #include <fevt/EyeTriggerData.h>
41 #include <fevt/Telescope.h>
42 #include <fevt/TelescopeSimData.h>
43 #include <fevt/TelescopeTriggerData.h>
44 #include <fevt/PixelSimData.h>
45 #include <fevt/Pixel.h>
46 #include <fevt/ChannelSimData.h>
47 #include <fevt/Channel.h>
48 
49 //++++++++++++++++++++++++++++++++++++++++
50 // Modified by D'Urso & Valore 02/03/06
51 #include <evt/LaserData.h>
52 //++++++++++++++++++++++++++++++++++++++++ end
53 
54 
55 // ROOT Headers
56 #include <TTree.h>
57 #include <TFile.h>
58 #include <TDirectory.h>
59 
60 
61 #include <io/FDasToOfflineEventConverter.h>
62 
63 #include "FdEventLogger.h"
64 
65 #include <AugerEvent.h>
66 
67 //-- FD-Event-Lib
68 #include <MiEvent.hh>
69 #include <MiEventHeader.hh>
70 #include <MiRun.hh>
71 #include <EyeEvent.hh>
72 
73 #include <FdNumbering.hh>
74 using namespace FdUtil;
75 
76 
77 using namespace FdEventLoggerGS;
78 using namespace std;
79 using namespace utl;
80 using namespace evt;
81 using namespace fwk;
82 using namespace det;
83 using namespace fevt;
84 using namespace fdet;
85 using namespace io;
86 
87 
88 
89 FdEventLogger::FdEventLogger() {}
90 FdEventLogger::~FdEventLogger() {}
91 
92 
94 
95  INFO ("FdEventLogger::Init()");
96 
97  CentralConfig* cc = CentralConfig::GetInstance();
98 
99  Branch topB = cc->GetTopBranch("FdEventLogger");
100  topB.GetChild("fileName").GetData(fFileName);
101 
102 
103  // Begin the output ROOT
104  TDirectory *const save = gDirectory;
105  fRootLogFile = new TFile(fFileName.c_str(), "RECREATE");
106  int comp = 1;
107  fRootLogFile->SetCompressionLevel(comp);
108 
109  // Init Root Tree
110  InitTree();
111  save->cd();
112 
113  return eSuccess;
114 
115 }// end of Init
116 
117 
118 
119 
120 VModule::ResultFlag FdEventLogger::Run(evt::Event& event){
121 
122  INFO("FdEventLogger::Run");
123 
124  if (!event.HasFEvent()) {
125  ERROR ("Event has no FEvent.");
126  return eContinueLoop;
127  }
128 
129  FEvent& fevent = event.GetFEvent();
130  int eventId = fevent.GetHeader().GetId();
131 
132  bool drumCalibMode = !event.HasSimShower();
133 
134  if (drumCalibMode) {
135  INFO ("Drum calibration mode!");
136  }
137 
138  primary = 0;
139  energy = 0;
140  XMax = 0;
141 
142  azimuth = 0;
143  zenith = 0;
144 
145  XG = 0;
146  YG = 0;
147 
148  laser_wavelength = 0;
149  laser_polarization_angle = 0;
150 
151  if (!drumCalibMode) {
152 
153  // Shower parameters for the Output Tree
154  evt::ShowerSimData& shower = event.GetSimShower();
155 
156  //++++++++++++++++++++++++++++++++++++++++
157  // Modified by D'Urso & Valore 02/03/06
158  bool LaserProcessing = shower.HasLaserData();
159  if (LaserProcessing == true) {
160  cout<<endl<<"...PROCESSING LASER EVENT..."<<endl;
161  laser_wavelength = shower.GetLaserData().GetLaserWavelength();
162  bool ispol = shower.GetLaserData().IsPolarized();
163  if (ispol == true) {
164  laser_polarization_angle = 0; //to be implemented
165  } else {
166  laser_polarization_angle = 0;
167  }
168  }
169  //++++++++++++++++++++++++++++++++++++++++ end
170 
171 
172 
173  //++++++++++++++++++++++++++++++++++++++++
174  // Modified by D'Urso & Valore 02/03/06
175  if (LaserProcessing==false) { //shower
176  energy = log10(shower.GetEnergy());
177  primary = shower.GetPrimaryParticle();
178  const evt::VGaisserHillasParameter& GHpointer = shower.GetGHParameters();
179  XMax = GHpointer.GetXMax() /(g/cm/cm);
180  } else { //laser
181  primary = 0;
182  energy = shower.GetEnergy()/utl::joule;
183  XMax = 0;
184  }
185  //++++++++++++++++++++++++++++++++++++++++ end
186 
187  if (shower.HasDirection()) {
188  const CoordinateSystemPtr localCS = shower.GetLocalCoordinateSystem();
189  zenith = (-shower.GetDirection()).GetTheta(localCS)/degree; // converted to degree
190  azimuth = (-shower.GetDirection()).GetPhi(localCS)/degree; // degree
191  }
192 
193  if (shower.HasPosition()) {
194 
195  const Point& core = shower.GetPosition();
196  const UTMPoint utm(core, ReferenceEllipsoid::eWGS84);
197  XG = utm.GetEasting();
198  YG = utm.GetNorthing();
199  }
200  }
201 
202 
203 
204  T2time = 0;
205  tslot = 1;
206 
207  fTimeOffset = 0;
208  fSLTWord[0] = 0;
209 
210 
211  for (FEvent::ConstEyeIterator iEye = fevent.EyesBegin(ComponentSelector::eHasData);
212  iEye != fevent.EyesEnd(ComponentSelector::eHasData);
213  ++iEye) {
214 
215  int eyeId = iEye->GetId();
216 
217  fIsT3 = false;
218  if (iEye->HasTriggerData()) {
219  fIsT3 = iEye->GetTriggerData().IsT3Accepted();
220  }
221 
222  for (fevt::Eye::ConstTelescopeIterator iTel = iEye->TelescopesBegin(ComponentSelector::eHasData);
223  iTel != iEye->TelescopesEnd(ComponentSelector::eHasData);
224  ++iTel) {
225 
226  int telId = iTel->GetId();
227  evtmir = 1000 * eventId + (eyeId - 1) * 6 + telId;
228 
229  if (!iTel->HasSimData()) {
230  INFO("Skipping telescope with no sim data");
231  continue;
232  }
233  const TelescopeSimData& telSim = iTel->GetSimData();
234  PixOn = telSim.GetNumberOfFltPixels();
235  PixOnFD = telSim.GetNumberOfFltPixelsFromShower();
236  PixOnBg = telSim.GetNumberOfFltPixelsFromBackground();
237  PixFromFD = telSim.GetNumberOfPixelsWithShowerPhotons();
238  PixRO = telSim.GetNumberOfReadOutPixels();
239  T2time = telSim.GetSltTriggerTime();
240  fTimeOffset = telSim.GetSltTimeShift();
241  tslot = 1;
242  fIsTLT = false;
243  if (iTel->HasTriggerData()) {
244  fIsTLT = iTel->GetTriggerData().IsTLTAccepted();
245  tslot = iTel->GetTriggerData().GetSLTData().size();
246  for(int i=0; i<tslot; i++) {
247  const int ncol = 20;
248  for(int col=0; col<ncol; col++) {
249  fSLTWord[i+ncol*col] = iTel->GetTriggerData().GetSLTData()[i].GetSLTDataWord(col+1);
250  }
251  }
252  }
253 
254 
255  // fill SLT tree
256  fT2trig->Fill();
257 
258 
259  for (fevt::Telescope::ConstPixelIterator iPixel = iTel->PixelsBegin(ComponentSelector::eHasData);
260  iPixel != iTel->PixelsEnd(ComponentSelector::eHasData);
261  ++iPixel) {
262 
263  int pixelId = iPixel->GetId();
264 
265  if (!iPixel->HasSimData()) {
266  INFO("Skip pixel with no sim data");
267  }
268  const PixelSimData& pixelSim = iPixel->GetSimData();
269 
270  // bool pixelTrigger = false;
271  // if (iPixel->HasTriggerData()) {
272  // pixelTrigger = true;
273  // } // unused
274 
275  pixel = pixelId;
276  thresh = pixelSim.GetThreshold();
277  summax = pixelSim.GetMaxBoxcarsum();
278  T1 = pixelSim.GetFltTime();
279  T1_delta = pixelSim.GetFltDuration();
280  mean = pixelSim.GetMeanBoxcarsum();
281  rms = pixelSim.GetRmsBoxcarsum();
282  FDSpix = pixelSim.HasPhotonTrace(fevt::FdConstants::eTotal);
283 
284 
285  //fill the pixel trigger tree
286  fTtrig->Fill();
287 
288  } // end loop over pixels
289 
290  } // end loop over Telescopes
291 
292  } // end loop over Eyes
293 
294  return eSuccess;
295 
296 }// end of RunDefault
297 
298 
299 
300 
301 VModule::ResultFlag FdEventLogger::Finish() {
302 
303  TDirectory* const save = gDirectory;
304 
305  INFO("Closing Analysis output Tree");
306  // Close Analysis output tree
307  fRootLogFile->cd();
308  fTtrig->Write();
309  fT2trig->Write();
310  fRootLogFile->Close();
311 
312  save->cd();
313 
314  return eSuccess;
315 
316 }// end of Finish
317 
318 
319 
320 void FdEventLogger::InitTree() {
321 
322  // Analysis output tree
323  TDirectory* const save = gDirectory;
324  fRootLogFile->cd();
325 
326  fTtrig = new TTree("pix","FirstLevelTrigger");
327  fTtrig->Branch("pixel",&pixel,"pixel/I");
328  fTtrig->Branch("evtmir",&evtmir,"evtmir/I");
329  fTtrig->Branch("thresh",&thresh,"thresh/I");
330  fTtrig->Branch("mean",&mean,"mean/F");
331  fTtrig->Branch("rms",&rms,"rms/F");
332  fTtrig->Branch("summax",&summax,"summax/I");
333  fTtrig->Branch("T1",&T1,"T1/I");
334  fTtrig->Branch("T1_delta",&T1_delta,"T1_delta/I");
335  fTtrig->Branch("FDSpix",&FDSpix,"FDSpix/b");
336 
337 
338  fT2trig = new TTree("SLT","SecondLevelTrigger");
339 
340  fT2trig->Branch("evtmir",&evtmir,"evtmir/I");
341  fT2trig->Branch("energy",&energy,"energy/D");
342  fT2trig->Branch("primary",&primary,"primary/I");
343  fT2trig->Branch("XMax",&XMax,"XMax/D");
344  fT2trig->Branch("azimuth",&azimuth,"azimuth/D");
345  fT2trig->Branch("zenith",&zenith,"zenith/D");
346  fT2trig->Branch("XG",&XG,"XG/D");
347  fT2trig->Branch("YG",&YG,"YG/D");
348  fT2trig->Branch("PixOn",&PixOn,"PixOn/I");
349  fT2trig->Branch("PixOnFD",&PixOnFD,"PixOnFD/I");
350  fT2trig->Branch("PixOnBg",&PixOnBg,"PixOnBg/I");
351  fT2trig->Branch("PixFromFD",&PixFromFD,"PixFromFD/I");
352  fT2trig->Branch("PixRO",&PixRO,"PixRO/I");
353  fT2trig->Branch("T2time",&T2time,"T2time/I");
354  fT2trig->Branch("T_offset",&fTimeOffset,"T_offset/I");
355  fT2trig->Branch("tslot",&tslot,"tslot/I");
356  fT2trig->Branch("SLTdataword",&fSLTWord,"SLTdataword[tslot]/i");
357  //++++++++++++++++++++++++++++++++++++++++
358  // Modified by D'Urso & Valore 02/03/06
359  fT2trig->Branch("LaserWavelength",&laser_wavelength,"LaserWavelength/D");
360  fT2trig->Branch("LaserPolarizationAngle",&laser_polarization_angle,"LaserPolarizationAngle/D");
361  //++++++++++++++++++++++++++++++++++++++++ end
362  fT2trig->Branch("isTLT", &fIsTLT, "isTLT/I");
363  fT2trig->Branch("isT3", &fIsT3, "isT3/I");
364 
365  save->cd();
366 }// end of InitTree
367 
368 
369 // Configure (x)emacs for this file ...
370 // Local Variables:
371 // mode:c++
372 // compile-command: "make -C .. -k"
373 // End:
bool HasDirection() const
Check initialization of shower geometry.
utl::CoordinateSystemPtr GetLocalCoordinateSystem() const
Get the Auger coordinate system associated to the shower core position.
int GetPrimaryParticle() const
Get the type of the shower primary particle.
Definition: ShowerSimData.h:84
LaserData & GetLaserData()
Get the laser data.
const double degree
Point object.
Definition: Point.h:32
bool HasLaserData() const
Check initialization of the LaserData.
double GetLaserWavelength() const
Definition: LaserData.h:32
Header & GetHeader()
Definition: FEvent.h:92
boost::filter_iterator< ComponentSelector, ConstAllEyeIterator > ConstEyeIterator
Definition: FEvent.h:56
const evt::VGaisserHillasParameter & GetGHParameters() const
Get the Gaisser-Hillas parameters of the shower.
bool HasFEvent() const
int GetId() const
Definition: FEvent/Header.h:21
int GetThreshold() const
Get the simulated trigger threshold of the pixel.
Definition: PixelSimData.h:132
Class to hold and convert a point in geodetic coordinates.
Definition: UTMPoint.h:40
int GetNumberOfFltPixelsFromBackground() const
EyeIterator EyesEnd(const ComponentSelector::Status status)
Definition: FEvent.h:66
#define INFO(message)
Macro for logging informational messages.
Definition: ErrorLogger.h:161
void Init()
Initialise the registry.
Branch GetChild(const std::string &childName) const
Get child of this Branch by child name.
Definition: Branch.cc:211
bool HasPhotonTrace(const FdConstants::LightSource source) const
Check that trace for source /par source is present.
Definition: PixelSimData.h:51
double GetNorthing() const
Get the northing.
Definition: UTMPoint.h:206
int GetNumberOfFltPixels() const
int GetSltTriggerTime() const
bool HasPosition() const
Check initialization of shower geometry.
Interface class to access Shower Simulated parameters.
Definition: ShowerSimData.h:49
int GetFltTime() const
Definition: PixelSimData.h:142
int GetMaxBoxcarsum() const
Definition: PixelSimData.h:139
boost::shared_ptr< const CoordinateTransformer > CoordinateSystemPtr
Shared pointer for coordinate systems.
Description of simulated data for one Telescope.
Class representing a document branch.
Definition: Branch.h:107
const utl::Vector & GetDirection() const
Get the direction of the shower axis. This is the true direction of shower movement.
const utl::Point & GetPosition() const
Get the position of the shower core.
EyeIterator EyesBegin(const ComponentSelector::Status status)
Definition: FEvent.h:58
boost::filter_iterator< ComponentSelector, ConstAllPixelIterator > ConstPixelIterator
double GetEasting() const
Get the easting.
Definition: UTMPoint.h:209
constexpr double g
Definition: AugerUnits.h:200
double GetEnergy() const
Get the energy of the shower primary particle.
Definition: ShowerSimData.h:89
int GetNumberOfPixelsWithShowerPhotons() const
void GetData(bool &b) const
Overloads of the GetData member template function.
Definition: Branch.cc:644
Top of Fluorescence Detector event hierarchy.
Definition: FEvent.h:33
constexpr double joule
Definition: AugerUnits.h:181
double GetMeanBoxcarsum() const
Definition: PixelSimData.h:140
int GetNumberOfReadOutPixels() const
ResultFlag
Flag returned by module methods to the RunController.
Definition: VModule.h:60
constexpr double cm
Definition: AugerUnits.h:117
bool IsPolarized() const
Definition: LaserData.h:38
Interface class for access to the Gaisser-Hillas parameters.
double GetRmsBoxcarsum() const
Definition: PixelSimData.h:141
int GetNumberOfFltPixelsFromShower() const
Main configuration utility.
Definition: CentralConfig.h:51
#define ERROR(message)
Macro for logging error messages.
Definition: ErrorLogger.h:165
boost::filter_iterator< ComponentSelector, ConstAllTelescopeIterator > ConstTelescopeIterator
Definition: FEvent/Eye.h:73
utl::Branch GetTopBranch(const std::string &id)
Get top branch for moduleConfigLink with given id (XML files)
Fluorescence Detector Pixel Simulated Data.
Definition: PixelSimData.h:28
int GetFltDuration() const
Definition: PixelSimData.h:143

, generated on Tue Sep 26 2023.