GammaWare  Head Version for release 0.9
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
BaseSpectrumPlayer.cpp
Go to the documentation of this file.
1 /***************************************************************************
2  * Copyright (C) 2004-2006 by Olivier Stezowski & Christian Finck *
3  * stezow(AT)ipnl.in2p3.fr, cfinck(AT)ires.in2p3.fr *
4  * *
5  * This program is free software; you can redistribute it and/or modify *
6  * it under the terms of the GNU General Public License as published by *
7  * the Free Software Foundation; either version 2 of the License, or *
8  * (at your option) any later version. *
9  * *
10  * This program is distributed in the hope that it will be useful, *
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
13  * GNU General Public License for more details. *
14  * *
15  * You should have received a copy of the GNU General Public License *
16  * along with this program; if not, write to the *
17  * Free Software Foundation, Inc., *
18  * 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. *
19  ***************************************************************************/
20 
23 #include <Riostream.h>
24 
25 
26 #ifndef ROOT_KeySymbols
27 #include "KeySymbols.h"
28 #endif
29 
30 #ifndef ROOT_TCanvas
31 #include "TCanvas.h"
32 #endif
33 
34 #ifndef ROOT_TClass
35 #include "TClass.h"
36 #endif
37 
38 #ifndef ROOT_TContextMenu
39 #include "TContextMenu.h"
40 #endif
41 
42 #ifndef ROOT_TH1
43 #include "TH1.h"
44 #endif
45 
46 #ifndef ROOT_TH2
47 #include "TH2.h"
48 #endif
49 
50 #ifndef ROOT_TFormula
51 #include "TFormula.h"
52 #endif
53 
54 #ifndef ROOT_TList
55 #include "TList.h"
56 #endif
57 
58 #ifndef ROOT_TMath
59 #include "TMath.h"
60 #endif
61 
62 #ifndef ROOT_TPad
63 #include "TPad.h"
64 #endif
65 
66 #ifndef ROOT_TROOT
67 #include "TROOT.h"
68 #endif
69 
70 #ifndef Gw_PeakCreator
71 #include "PeakCreator.h"
72 #endif
73 
74 #ifndef Gw_BaseSpectrumPlayer
75 #include "BaseSpectrumPlayer.h"
76 #endif
77 
78 
79 using namespace Gw;
80 
81 //
82 
83 //__________________________________________________________
85  : TNamed("SpectrumPlayer","Player for spectra"),
86  fCreator(0x0),
87  fPeakList(0x0),
88  fInnerPeakList(0x0),
89  fLastX(0),
90  fLastY(0),
91  fIsHelpsPrintActive(true),
92  fContextMenu(0x0),
93  fLog("BaseSpectrumPlayer")
94 {
95  fCreator = new PeakCreator();
96 
97  // default constructor
98  fInnerPeakList = new TList(); fInnerPeakList->SetOwner(true);
99  fPeakList = fInnerPeakList;
100 }
101 
102 //__________________________________________________________
104 {
105  // default destructor
106  fPeakList = 0x0;
107 
108  if (fCreator)
109  { delete fCreator; fCreator = 0x0; }
110 
111  fInnerPeakList->Delete();
112  if (fInnerPeakList)
113  { delete fInnerPeakList; fInnerPeakList = 0x0; }
114 
115  Disconnect();
116 }
117 
118 
119 //__________________________________________________________
121 {
122  if ( h->GetDimension() != 1) {
123  fLog.SetProcessMethod("IsInRange(TH1* , BasePeak* )");
124  fLog << error << "2D peak method not implemented yet" << dolog;
125  return false;
126  }
127 
128  Int_t firstBin = h->GetXaxis()->GetFirst();
129  Int_t lastBin = h->GetXaxis()->GetLast();
130 
131  Float_t low = h->GetBinCenter(firstBin);
132  Float_t high = h->GetBinCenter(lastBin);
133 
134  if (peak->GetPosition() - 2*peak->GetFWHM() > low && peak->GetPosition() + 2*peak->GetFWHM() < high)
135  return true;
136  else {
137  fLog.SetProcessMethod("IsInRange(TH1* , BasePeak* )");
138  fLog << info << Form("Peak %s not fitted cos out of zoom", peak->GetName()) << dolog;
139  return false;
140  }
141 
142  return false;
143 }
144 
145 //__________________________________________________________
146 Bool_t BaseSpectrumPlayer::Connect(TCanvas *canvas)
147 {
148  TCanvas *localCanvas = canvas;
149 
150  // in this case connect the current canvas
151  if ( canvas == 0x0 ) {
152  if ( TVirtualPad::Pad() )
153  localCanvas = TVirtualPad::Pad()->GetCanvas();
154  else {
155  fLog.SetProcessMethod("Connect(TCanvas* )");
156  fLog << error << "gPad not existing" << dolog;
157  }
158  }
159  if ( localCanvas) {
160 
161  // to fit graphically
162  localCanvas->Connect("ProcessedEvent(Int_t, Int_t, Int_t, TObject*)",
163  "Gw::BaseSpectrumPlayer",
164  this,
165  "HandleMovement(Int_t, Int_t, Int_t, TObject*)");
166 
167  fContextMenu = localCanvas->GetContextMenu();
168 
169  // add a peak creator to the canvas if it does not exists
170 
171  if ( ! PeakCreator::IsPeakCreator(localCanvas) ) {
172  // change the characteristics ?
173  // PeakCreator *p = PeakCreator::AddPeakCreator(localCanvas);
175  fCreator->Connect();
176  }
177  else
178  {
179  fCreator->Connect();
180  }
181  fCreator->Collect();
182  }
183  return localCanvas != 0x0;
184 }
185 
186  //__________________________________________________________
187 Bool_t BaseSpectrumPlayer::Disconnect(TCanvas *canvas)
188 {
189  fLog.SetProcessMethod("Disconnect(TCanvas* )");
190 
191 
192  TCanvas *localCanvas = canvas;
193 
194  // in this case connect the current canvas
195  if ( canvas == 0x0 ) {
196  if ( TVirtualPad::Pad() )
197  localCanvas = TVirtualPad::Pad()->GetCanvas();
198  else {
199  fLog << error << "gPad not existing" << dolog;
200  }
201  }
202  if ( localCanvas ) {
203  if (localCanvas->Disconnect("ProcessedEvent(Int_t, Int_t, Int_t, TObject*)",
204  this,
205  "HandleMovement(Int_t, Int_t, Int_t, TObject*)"))
206  fLog << info << "Disconnected from Canvas" << dolog;
207  else
208  fLog << error << "Cannot disconnected from Canvas" << dolog;
209 
210  }
211  fContextMenu = 0x0;
212 
213  PeakCreator *PC = PeakCreator::IsPeakCreator(localCanvas);
214  if(PC) PC->Disconnect();
215 
216  fLog << dolog;
217  return localCanvas != 0x0;
218 }
219 
220 //__________________________________________________________
221 void BaseSpectrumPlayer::HandleMovement(Int_t eventType, Int_t eventX, Int_t eventY, TObject* /*select*/)
222 {
223 // // handle cursor mouvement
224 // TVirtualPad *pad = gROOT->GetSelectedPad(), *cpad = TVirtualPad::Pad();
225 
226 // if ( pad != cpad )
227 // pad->cd();
228 
229  if( eventType == kKeyPress) {
230 
231  Int_t keysym = eventY;
232 
233  if((EKeySym)eventY == kKey_a && (EKeySym)eventX == kKey_a && (EKeySym)fLastY == kKey_f && (EKeySym)fLastX == kKey_f)
234  PopupFitMenu();
235 
236  if((EKeySym)eventY == kKey_m && (EKeySym)eventX == kKey_m && (((EKeySym)fLastY == kKey_s && (EKeySym)fLastX == kKey_s) || ((EKeySym)fLastY == kKey_h && (EKeySym)fLastX == kKey_h)))
237  {
238  if(TVirtualPad::Pad() == 0x0) return;
239  TList *list = TVirtualPad::Pad()->GetListOfPrimitives();
240  if(list==0x0) return;
241 
242  for(int i=0 ; i<list->GetSize() ; i++)
243  {
244  TObject *o = list->At(i);
245  if(o->InheritsFrom("Gw::Peak1D"))
246  {
247  Peak1D *peak = (Peak1D*)o;
248  if(peak->IsDrawAs(Peak1D::kGate)) continue;
249  else
250  {
251  if((EKeySym)fLastY == kKey_s) peak->SetDrawOption("ml");
252  else peak->SetDrawOption("!m!l");
253  }
254  }
255  }
256 
257  gPad->Modified();
258  gPad->Update();
259  }
260 
261  if ( (EKeySym)keysym == kKey_h && fIsHelpsPrintActive) {
262  fLog << info;
263  fLog << "\tf + a: pop up fit panel menu" << nline;
264  fLog << dolog;
265  }
266  }
267 
268  fLastX = eventX; fLastY = eventY;
269 
270 // cpad->cd();
271 }
272 
273 //__________________________________________________________
275 {
276  TMethod* m = Class()->GetMethod("FitAll", "");
277 
278  if (m != 0x0 && fContextMenu != 0x0)
279  fContextMenu->Action(this, m);
280  else {
281  fLog.SetProcessMethod("PopupFitMenu()");
282  fLog << error << "Cannot Popup FitAll dialog box" << dolog;
283  }
284 }
285 
286 
287 //__________________________________________________________
288 void BaseSpectrumPlayer::FitAll(const char* nameFunc, Option_t* optFit, Option_t* optBkg)
289 {
290  // fit all peak in list with a common function and background
291  if ( CollectPeaks() == 0 )
292  return;
293 
294  BasePeak* peak = 0x0;
295  BasePeak* peak2 = 0x0;
296 
297  TString FuncName = nameFunc;
298 // FuncName.ToLower();
299  TString BackOpt = optBkg;
300 // BackOpt.ToLower();
301 
302  if(FuncName == "DTGaus" && BackOpt.Contains("Step"))
303  {
304  FuncName = "DTGaus_Step";
305  BackOpt = "";
306  }
307 
308  TList peaks; peaks.SetOwner(false);
309 
310  // sort peak in list
311  SortPeakList();
312 
313  TString opt(optFit);
314  // rename if name of peaks are the same
315  if ( RenamePeak(GetName()) ) {
316  fLog.SetProcessMethod("FitAll(const char*, Option_t*, Option_t* )");
317  fLog << warning << "Rename peaks cos they have the same name" << dolog;
318  }
319 
320  TH1 *h = PadManager::GetHisto();
321 
322  if(h == 0x0)
323  return;
324 
325  if ( h->GetDimension() >= 2) {
326  fLog.SetProcessMethod("FitAll(const char*, Option_t*, Option_t* )");
327  fLog << error << "2D peak fitting not implemented yet" << dolog;
328  return;
329  }
330 
331  Bool_t IsUsed[GetPeakList()->GetEntries()];
332  for (Int_t i = 0; i < GetPeakList()->GetEntries(); ++i)
333  IsUsed[i] = false;
334 
335  for (Int_t i = 0; i < GetPeakList()->GetEntries(); ++i) {
336  if (!IsUsed[i]) IsUsed[i] = true;
337  else continue;
338 
339  peak = static_cast<BasePeak*> (GetPeakList()->At(i));
340  if (opt.Contains("R") && !IsInRange(h, peak)) continue;
341 
342  for (Int_t j = i+1; j < GetPeakList()->GetEntries(); ++j) {
343  peak2 = static_cast<BasePeak*> (GetPeakList()->At(j));
344  if (opt.Contains("R") && !IsInRange(h, peak2)) continue;
345 
346  if (peak->IsPointInBkg(peak2->GetPosition()) < 1) {
347  peak2->SetSignalFunction(FuncName,h);
348  if (!IsUsed[j]) {
349  peaks.Add(peak2);
350  IsUsed[j] = true;
351  }
352  }
353  }
354 
355  if (peaks.GetEntries() != 0)
356  {
357  if(FuncName == "DTGaus_Step")
358  {
359  FuncName = "DTGaus";
360  BackOpt = "lin";
361  }
362 
363  peak->SetSignalFunction(FuncName,h);
364  for(int j=0 ; j<peaks.GetEntries() ; j++)
365  ((BasePeak*)peaks.At(j))->SetSignalFunction(FuncName,h);
366 
367  peak->FitCombined(h, &peaks, optFit, BackOpt);
368  peaks.Clear("nodelete");
369  } else {
370 
371  peak->SetSignalFunction(FuncName,h);
372  peak->Fit(h,optFit, BackOpt);
373  }
374  }
375 
376  TIter next(&peaks); TObject *obj;
377  while ( (obj = next()) ) {
378  peaks.Remove(obj);
379  }
380 
381  gPad->Modified();
382  gPad->Update();
383 }
384 
385 
386 //__________________________________________________________
388 {
389  // set bkg from current histo in pad
390  TH1 *h = PadManager::GetHisto(), *result = 0x0;
391  if ( h )
392  result = Background(h,opt);
393  if ( result ) {
394  result->Draw("same");
395  } // TODO Wheel Color/Style
396 
397  return result;
398 }
399 
400 //__________________________________________________________
401 void BaseSpectrumPlayer::DoBackground(TH1 *histo, Option_t* opt)
402 {
403  // set bkg from histo
404  TH1 *bg = Background(histo,opt);
405  if ( bg ) {
406  histo->Add(bg,-1.0);
407  delete bg;
408  }
409 }
410 
411 //__________________________________________________________
412 void BaseSpectrumPlayer::SetPeakList(TSeqCollection *col)
413 {
414  // set peak list from outside
415  col ? fPeakList = col : fPeakList = dynamic_cast<TSeqCollection*>(fInnerPeakList);
416 }
417 
418 //__________________________________________________________
419 Bool_t BaseSpectrumPlayer::SetParameter(const char* name, const TObject *value)
420 {
421  fLog.SetProcessMethod("SetParameter(const char*, TObject* )");
422 
423  // set parameters from list
424  Bool_t ok = false;
425 
426  if ( value == 0x0 )
427  return ok;
428 
429  // set search parameter
430  if (strncmp(name, "DefaultPeak1D", strlen(name)) == 0 && value->InheritsFrom("Gw::BasePeak") ) {
431 /* if ( fDefaultPeak1D )
432  { delete fDefaultPeak1D; fDefaultPeak1D = 0x0; }
433  fDefaultPeak1D =
434  dynamic_cast<BasePeak *>(value->Clone());
435  ok = true; */
436  }
437  else {
438  fLog << error << "Unkown parameter " << name << nline;
439  }
440 
441  fLog << dolog; return ok;
442 }
443 
444 //__________________________________________________________
446 {
447  PeakCreator* peakCreator = PeakCreator::IsPeakCreator(TVirtualPad::Pad());
448  if (peakCreator) {
449  peakCreator->SetDefaultPeakFWHM(name);
450  } else {
451  fLog.SetProcessMethod("SetDefaultPeakFWHM(const char* )");
452  fLog << warning << "Peak Creator not existing in the current Pad" << dolog;
453  }
454 
455 }
456 
457 //__________________________________________________________
458 Int_t BaseSpectrumPlayer::FindPeaks(Option_t* opt)
459 {
460  // find peak method
461  TH1* h = PadManager::GetHisto();
462  if ( h == 0x0 ) {
463  fLog.SetProcessMethod("FindPeaks(Option_t* )");
464  fLog << error << "No spectrum displayed in current pad" << dolog;
465  return 0;
466  }
467  else return FindPeaks(h,opt);
468 }
469 
470 //__________________________________________________________
471 Bool_t BaseSpectrumPlayer::RenamePeak(const Char_t* baseName, Bool_t force)
472 {
473  // rename and re-ordered peak with respect to position
474  BasePeak* peak = 0x0;
475  BasePeak* peak2 = 0x0;
476  Int_t k = 0;
477 
478  Bool_t rename = false;
479  for (Int_t i = 0; i < GetPeakList()->GetEntries(); ++i) {
480  peak = static_cast<BasePeak*> (GetPeakList()->At(i));
481  TString name(peak->GetName());
482  for (Int_t j = i+1; j < GetPeakList()->GetEntries(); ++j) {
483  peak2 = static_cast<BasePeak*> (GetPeakList()->At(j));
484 
485  if (name.CompareTo(peak2->GetName()) == 0)
486  rename = true;
487  }
488  }
489 
490  if (!force && !rename) return false;
491 
492  TIter next(fPeakList);
493  while ( (peak = (BasePeak*) next()) ) { // iterator skips empty slots
494  peak->SetName( Form("%s%d", baseName, k++) );
495  }
496  return true;
497 }
498 
499 //__________________________________________________________
501 {
502  // add peak to list
503  if ( !fPeakList->FindObject(peak) )
504  fPeakList->Add(peak);
505 }
506 
507 //__________________________________________________________
508 void BaseSpectrumPlayer::SortPeakList(const Char_t* parName, Bool_t sortDes)
509 {
510  // sort peak respect to given parameter
511 
512  BasePeak* peak = 0x0;
513  Int_t i = 0;
514  Int_t n = fPeakList->GetEntries();
515 
516  Int_t* index = new Int_t[n];
517  Double_t* par = new Double_t[n];
518 
519  TList* tmp = new TList();
520  tmp->SetOwner(false);
521 
522  // save parameters
523  TIter next(fPeakList);
524  while ( (peak = (Gw::BasePeak*)next()) ) {
525  if (strncmp(parName, "Position", n) == 0)
526  par[i++] = peak->GetPosition();
527  if (strncmp(parName, "FWHM", n) == 0)
528  par[i++] = peak->GetFWHM();
529  if (strncmp(parName, "Height", n) == 0)
530  par[i++] = peak->GetHeight();
531  tmp->AddLast(peak);
532  }
533 
534  // sort in index
535  TMath::Sort(n, par, index, sortDes);
536 
537  // just reset counter, no delete
538  fPeakList->Clear("nodelete");
539 
540  // add again in fInnerPeakList
541  for (i = 0; i < n; ++i ) {
542  peak = static_cast<BasePeak*>(tmp->At(index[i]));
543  fPeakList->AddLast(peak);
544  }
545 
546  tmp->Clear("nodelete");
547  delete tmp;
548  delete[] index;
549  delete[] par;
550 }
551 
552 //__________________________________________________________
554 {
555  // re-draw peaks onto current pad
556  fPeakList->AppendPad();
557 }
558 
559 //__________________________________________________________
561 {
562  // create peak list from pad
563  if ( !TVirtualPad::Pad() ) {
564  fLog.SetProcessMethod("CollectPeaks(Option_t* )");
565  fLog << error << "Pad not existing" << dolog;
566  return 0;
567  }
568 
569  TString opt = o;
570  if ( !opt.Contains("+") )
571  fPeakList->Clear("nodelete");
572 
573  TObject *obj; Int_t i = 0;
574 
575  TIter next(TVirtualPad::Pad()->GetListOfPrimitives());
576  while ( (obj = next()) ) { // iterator skips empty slots
577  if ( obj->InheritsFrom("Gw::BasePeak") ) {
578  if (!fPeakList->FindObject(obj)) {
579 // obj->SetBit(TObject::kCanDelete, false);
580 // obj->SetBit(TObject::kMustCleanup, false);
581  fPeakList->Add(obj); i++;
582  }
583  } }
584 
585  return i;
586 }
587 
588 //__________________________________________________________
589 void BaseSpectrumPlayer::Print(Option_t* opt) const
590 {
591  // print peak info
592  fPeakList->Print(opt);
593 }
594 
596 
void SetDefaultPeakFWHM(const char *fname="PCG_FWHM")
void SetDefaultPeakFWHM(const char *="PCG_FWHM")
to change the definition of the formula for the width
virtual Double_t GetPosition(Option_t *axis="X") const =0
Get position of peak.
header file for PeakCreator.cpp
virtual TF1 * SetSignalFunction(const char *nameFunc="gaus", TH1 *h=0x0)=0
Set pre-defined function to fit the signal.
virtual void SortPeakList(const Char_t *parName="Position", Bool_t sortDes=false)
sort peak list
virtual void Fit(TH1 *h, Option_t *optFit="RN", Option_t *optBkg="lin")=0
Fit peak with background.
LogMessage & error(LogMessage &)
LogMessage & warning(LogMessage &)
virtual void AddCreatorToList()
Definition: PeakCreator.h:98
virtual void DoBackground(TH1 *histo, Option_t *opt="")
Apply background substraction for that histogram.
BaseSpectrumPlayer to work on spectra.
virtual Bool_t RenamePeak(const Char_t *baseName="Peak", Bool_t force=false)
Sort all peaks and rename them.
LogMessage fLog
log message
BaseSpectrumPlayer()
default constructor
Bool_t Disconnect(TCanvas *c=0x0)
Connect the Canvas to this to collect events.
bool IsDrawAs(EDrawAs mode)
To test the draw mod ex: IsDrawAs(kGate)
Definition: Peak1D.h:366
LogMessage & nline(LogMessage &)
UInt_t value[MaxValue]
Definition: ReadDaqAlone.C:29
TH1F * histo[MaxValue]
Definition: ReadDaqAlone.C:31
virtual Double_t GetHeight() const
Get height of peak.
Definition: BasePeak.h:66
virtual void Connect()
header file for BaseSpectrumPlayer.cpp
virtual Short_t IsPointInBkg(Double_t x, Double_t y=0)=0
to determine if a point is in bg. 0 likely in peak, 1 in bg, 2 outside bg
virtual Int_t FindPeaks(Option_t *opt="")
Find peaks from the histo in the current pad and store them in the current collection.
virtual void PopupFitMenu()
Popup AddLink menu.
LogMessage & info(LogMessage &)
manipulator to modify the LogMessage
virtual void FitCombined(TH1 *h, TList *peakList, Option_t *optFit="RN", Option_t *optBkg="lin")=0
Fit multi-defined peak with a common background.
PeakCreator * fCreator
default PeakCreator
TContextMenu * fContextMenu
context menu in canvas
void Collect(Bool_t do_collect=true)
if true, peak are collected on key board actions (type h for help)
LogMessage & dolog(LogMessage &)
void HandleMovement(Int_t eventType, Int_t eventX, Int_t eventY, TObject *)
Handle Movement.
Bool_t IsInRange(TH1 *h, BasePeak *peak)
check if peak is in range of drawn spectrum
ClassImp(BaseNucleus)
virtual void SetPeakList(TSeqCollection *col=0x0)
set collection
virtual Double_t GetFWHM(Option_t *axis="X") const =0
Get FWHM of peak.
virtual void Print(Option_t *opt="") const
Print out informations concerning the parameters of that player.
virtual TSeqCollection * GetPeakList() const
Return the current collection of peaks.
virtual void AddPeak(BasePeak *peak)
Add a peak to the current collection (at the end)
virtual Int_t CollectPeaks(Option_t *opt="")
Collect the peaks from the current pad.
A BasePeak is defined by a height, intensity and a dimension of the peak.
Definition: BasePeak.h:19
static PeakCreator * IsPeakCreator(TVirtualPad *pad=0x0)
Check whether or not a PeakCreator is in this canvas.
virtual void Disconnect()
virtual void FitAll(const char *nameFunc="DTGaus", Option_t *optFit="RN", Option_t *optBkg="Step")
Fit all peaks using the spectrum in the current pad.
A graphical interface for placing schematic peak onto a 1D histogram with a given position...
Definition: Peak1D.h:79
virtual void ShowPeakList() const
Show the list of peaks on the current pad.
virtual void SetProcessMethod(const char *)
To set the current method.
static TH1 * GetHisto(TVirtualPad *pad=0x0, Option_t *op="")
look for an histogram into the pad
Definition: PadManager.cpp:67
virtual TH1 * Background(const TH1 *, Option_t *)
Compute the background for that histogram.
virtual Bool_t SetParameter(const char *, Double_t)
To change the parameters for that algorithm.
Bool_t Connect(TCanvas *c=0x0)
Connect the Canvas to this to collect events.
virtual void SetDrawOption(Option_t *option="")
Definition: Peak1D.cpp:1796