NCFS-Pack
A generic (astro)particle physics analysis framework
 
Loading...
Searching...
No Matches
RnoStation.cxx
Go to the documentation of this file.
1
17
19
41
42#include "RnoStation.h"
43#include "Riostream.h"
44
45ClassImp(RnoStation); // Class implementation to enable ROOT I/O
46
48RnoStation::RnoStation(const char* name,const char* title) : NcDetectorUnit(name,title)
49{
55
56 fCanvas=0;
57}
58
60{
66
67 if (fCanvas)
68 {
69 if (gROOT->GetListOfCanvases()->FindObject(ClassName())) delete fCanvas;
70 fCanvas=0;
71 }
72}
73
84
85TGraph* RnoStation::DisplaySampling(Int_t ich,Int_t j)
86{
98
99 if (fCanvas)
100 {
101 if (gROOT->GetListOfCanvases()->FindObject(ClassName())) delete fCanvas;
102 fCanvas=0;
103 }
104
105 if (ich<0 || j<1) return 0;
106
107 fCanvas=new TCanvas();
108 TVirtualPad* pad=fCanvas->cd();
109 gROOT->SetSelectedPad(pad);
110
111 TString devname="Ch";
112 devname+=ich;
113
114 Int_t id=GetUniqueID(); // The station ID
115 TString staname="Station";
116 staname+=id;
117
118 TString title="";
119
120 title=staname;
121 title+=" ";
122 title+=devname;
123 fCanvas->SetName(ClassName());
124 fCanvas->SetTitle(title);
125 fCanvas->SetGrid();
126
127 GetSamplingGraph(ich,j).DrawClone("AL");
128
129 // Obtain the pointer to the displayed graph
130 TGraph* gr=(TGraph*)fCanvas->FindObject("NcSample");
131
132 return gr;
133}
134
136{
148
149 if (fCanvas)
150 {
151 if (gROOT->GetListOfCanvases()->FindObject(ClassName())) delete fCanvas;
152 fCanvas=0;
153 }
154
155 if (j<1) return 0;
156
157 fCanvas=new TCanvas();
158
159 Int_t id=GetUniqueID(); // The station ID
160 TString staname="Station";
161 staname+=id;
162
163 TString title="";
164 TString devname="";
165
166 title=staname;
167 fCanvas->SetName(ClassName());
168 fCanvas->SetTitle(title);
169 fCanvas->Divide(4,6);
170
171 TVirtualPad* pad=0;
172 for (Int_t jch=0; jch<24; jch++)
173 {
174 pad=fCanvas->cd(jch+1);
175
176 if (!pad) continue;
177
178 gROOT->SetSelectedPad(pad);
179 pad->SetGrid();
180
181 GetSamplingGraph(jch,j).DrawClone("AL");
182 } // End of loop over the channels
183
184 return fCanvas;
185}
186
187TGraph RnoStation::GetSamplingGraph(Int_t ich,Int_t j)
188{
197
198 TGraph gr;
199
200 if (ich<0 || j<1) return gr;
201
202 TString devname="Ch";
203 devname+=ich;
204
205 Int_t id=GetUniqueID(); // The station ID
206 TString staname="Station";
207 staname+=id;
208
209 NcDevice* daq=GetDevice("DAQ",kFALSE);
210 Double_t fsample=0;
211 if (daq) fsample=daq->GetSignal("Sampling-rate");
212
213 TString rate="";
214 rate.Form(" (DAQ: %-.3g Samples/sec)",fsample);
215
216 TString title=staname;
217 title+=" ";
218 title+=devname;
219 title+=rate;
220
221 RnoGANT* ant=(RnoGANT*)GetDevice(devname,kTRUE);
222 NcSample* sx=ant->GetSample(1);
223
224 if (!sx) return gr;
225
226 Int_t ndim=sx->GetDimension();
227 if (j>ndim) return gr;
228
229 gr=sx->GetGraph(j);
230
231 TString xname="Sample";
232 TString yname=sx->GetVariableName(j);
233
234 gr.SetTitle(title);
235 gr.GetXaxis()->SetTitle(xname);
236 gr.GetYaxis()->SetTitle(yname);
237
238 return gr;
239}
240
241TH1F RnoStation::GetSamplingDFT(Int_t ich,TString sel,Int_t j)
242{
269
270 TH1F his;
271
272 TGraph gr=GetSamplingGraph(ich,j);
273
274 NcDevice* daq=GetDevice("DAQ",kFALSE);
275 Double_t fsample=0;
276 if (daq) fsample=daq->GetSignal("Sampling-rate");
277
278 NcDSP q;
279 q.Load(&gr,fsample);
280 q.Fourier("R2C",&his,sel);
281
282 TString title=gr.GetTitle();
283
284 his.SetTitle(title);
285 his.SetStats(0);
286
287 return his;
288}
289
290TObject* RnoStation::Clone(const char* name) const
291{
301
302 RnoStation* q=new RnoStation(*this);
303 if (name)
304 {
305 if (strlen(name)) q->SetName(name);
306 }
307 return q;
308}
309
ClassImp(RnoStation)
Various Digital Signal Processing (DSP) operations for (sequential) data samples.
Definition NcDSP.h:21
void Fourier(TString mode, TH1 *hist=0, TString sel="none")
Definition NcDSP.cxx:1092
void Load(Int_t n, Double_t *re, Double_t *im=0, Float_t f=-1)
Definition NcDSP.cxx:395
NcDetectorUnit(const char *name="", const char *title="")
NcDevice(const char *name="", const char *title="")
Definition NcDevice.cxx:109
Sampling and statistics tools for various multi-dimensional data samples.
Definition NcSample.h:28
TString GetVariableName(Int_t i) const
Int_t GetDimension() const
TGraph GetGraph(Int_t i, TF1 *f=0)
NcDevice * GetDevice() const
NcSample * GetSample(Int_t j=1) const
virtual Float_t GetSignal(Int_t j=1, Int_t mode=0) const
Definition NcSignal.cxx:651
Signal (Hit) handling of an RNO-G Generic Antenna (GANT).
Definition RnoGANT.h:12
Handling of RNO-G event data.
Definition RnoStation.h:19
TH1F GetSamplingDFT(Int_t ich, TString sel="AMP Hz", Int_t j=1)
TCanvas * fCanvas
! Pointer to the temp. canvas for displays
Definition RnoStation.h:31
TGraph GetSamplingGraph(Int_t ich, Int_t j=1)
TCanvas * DisplaySamplings(Int_t j=1)
virtual TObject * Clone(const char *name="") const
virtual ~RnoStation()
RnoStation(const char *name="", const char *title="")
TGraph * DisplaySampling(Int_t ich, Int_t j=1)