GammaWare  Head Version for release 0.9
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GANIL/TrackedWatchers.C
Go to the documentation of this file.
1 /***************************************************************************
2  * Copyright (C) 2010 by Olivier Stezowski *
3  * stezow(AT)ipnl.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 
21 // Watcher definition
22 #include "TrackedWatchers.h"
23 
24 // #include "GenericFrame.h"
25 #include "AgataKeyFactory.h"
26 
27 #include "MetaWatchers.h"
28 #include "TVector3.h"
29 #include "TGeoMatrix.h"
31 
32 #include "TPaveStats.h"
33 #include "TLegend.h"
34 
35 //#define APPLY_CORRECTION_CrystalPositionLookUpTable
36 
37 using namespace ADF;
38 
42 
44 {
45  Bool_t ok = Watcher::GetFromTrigger(trigger, "data:tracked", fFrame);
46  return ok;
47 }
48 
50 
51 void TrackedWatcher::SetAgataOffset(Double_t offset)
52 {
53  fZOffset = offset;
54 }
55 
57 
59 {
60 #ifdef APPLY_CORRECTION_CrystalPositionLookUpTable
61  // ~30 degre in CrystalPositionLookupTable
62 
63  // cout << " inside Tracked watchers " << endl;
64 
65  Double_t e = hit->GetE(), beta = vertex->GetBeta(), gamma = vertex->GetGamma();
66 
67  Double_t px, py, pz, dx, dy, dz;
68  vertex->GetPosition (px, py, pz);
69  vertex->GetDirection(dx, dy, dz);
70 
71  Double_t hrx, hry, hrz; // rotated hitx,y,z
72  Double_t phi = -28*3.141592654/180; // rotation around z axis of ?? degrees
73  // Double_t phiy = 1.*3.141592654/180; // rotation around z axis of ?? degrees
74 
75  // cout << cos(phix) << " " << sin(phix) << endl;
76 
77 
78  // hrx = hit->GetX();
79  // hry = hit->GetY()*cos(phix) - hit->GetZ()*sin(phix);
80  // hrz = hit->GetY()*sin(phix) + hit->GetZ()*cos(phix);
81  hrx = hit->GetX()*cos(phi) - hit->GetY()*sin(phi);
82  hry = hit->GetX()*sin(phi) + hit->GetY()*cos(phi);
83  hrz = hit->GetZ();
84  Double_t xx = hit->GetX() - px;
85  Double_t yy = hit->GetY() - py;
86  Double_t zz = hit->GetZ() + fZOffset - pz;
87 
88  // Double_t x = hit->GetX() - px;
89  // Double_t y = hit->GetY() - py;
90  // Double_t z = hit->GetZ() + fZOffset - pz;
91  Double_t x = hrx - px;
92  Double_t y = hry - py;
93  Double_t z = hrz + fZOffset - pz;
94  Double_t dd = x*x + y*y + z*z;
95  Double_t ddd = xx*xx + yy*yy + zz*zz;
96 
97  Double_t enr = 0.;
98 
99  if( !(dd == 0.0)&&beta>0. ) {
100  Double_t cosTheta = (x*dx + y*dy + z*dz)/::sqrt(dd);
101  Double_t cosThetanr = (xx*dx + yy*dy + zz*dz)/::sqrt(ddd); //cos theta if not rotating by phix, for tests
102  Double_t edop=e;
103 
104  e = VertexInterface::DopplerCorrection(e,beta,gamma,cosTheta);
105  // enr = VertexInterface::DopplerCorrection(edop,beta,gamma,cosThetanr);
106 
107  if( abs(e - 1345)<10. && 0==-1 ) {
108  cout << " e before =" << edop << "; beta =" << beta << endl;
109  cout << " rotation: " << e << " " << cosTheta << "; no rotation: " << enr << " " << cosThetanr << " " << endl;
110  // cout << " with rotation: " << e << " " << cosTheta << " " << endl;
111 
112  cout << "check " << hit->GetX() << " --> " << hrx << " " << endl;
113  cout << "check " << hit->GetY() << " --> " << hry << " " << endl;
114  cout << "check " << hit->GetZ() << " --> " << hrz << " " << endl;
115  cout << "check " << hit->GetT() << endl;
116  // cout << "check " << xx << " --> " << x << " " << endl;
117  // cout << "check " << yy << " --> " << y << " " << endl;
118  // cout << "check " << zz << " --> " << z << " " << endl;
119  cout << " " << endl;
120  }
121 
122  }
123  else e = -1.0;
124 
125  return e;
126 
127 #else
128  Double_t e = hit->GetE(), beta = vertex->GetBeta(), gamma = vertex->GetGamma();
129 
130  Double_t px, py, pz, dx, dy, dz;
131  vertex->GetPosition (px, py, pz);
132  vertex->GetDirection(dx, dy, dz);
133 
134  Double_t x = hit->GetX() - px;
135  Double_t y = hit->GetY() - py;
136  Double_t z = hit->GetZ() + fZOffset - pz;
137 
138  Double_t dd = x*x + y*y + z*z;
139  if( !(dd == 0.0) ) {
140  Double_t cosTheta = (x*dx + y*dy + z*dz)/::sqrt(dd);
141  e = VertexInterface::DopplerCorrection(e,beta,gamma,cosTheta);
142  }
143  else e = -1.0;
144 
145  return e;
146 #endif
147 }
148 
149 
153 
154 Coinc2D::Coinc2D(const char *name, const char *title, TDirectory *sp_dir, TDirectory *tag_dir):
155  TrackedWatcher(name,title,sp_dir,tag_dir),
156  fVertexBuilder(0x0),
157  fTrackedEnergy(0x0),
158  fTrackedEnergyDopler(0x0),
159  fFold(0x0),
160  fGxG(0x0),
161  fGxG_DC(0x0),
162  fPos3D(0x0)
163 {
164  // Vertex used for doppler correction, last registered
165  fVertexBuilder = VertexBuilder::theCurrentVertexBuilder();
166 
167  // Now spectra
168  fTrackedEnergy = MakeTH1<TH1F>("TrackedEnergy","Tracked Energy",6000,0,6000);
169 
170  fTrackedEnergyDopler = MakeTH1<TH1F>("TrackedEnergyDopler","Tracked Energy with Doppler correction",6000,0,6000);
171  TagOn(fTrackedEnergyDopler);
172 
173  fFold = MakeTH1<TH1F>("Fold","Fold distributions",10,0,10);
174  TagOn(fFold);
175 
176  fGxG = MakeTH2<TH2F>("GxG","Gamma Gamma matrix",2000,0,1000,2000,0,1000);
177  fGxG_DC = MakeTH2<TH2F>("GxG_DC","Gamma Gamma matrix with Doppler Correction",2000,0,1000,2000,0,1000);
178 
179  fPos3D = MakeTH3<TH3F>("Pos3D","Positions of Tracked Gamma-rays ;x;y;z",200,-200,200,200,-200,200,200,-400,0);
180 }
181 
183 
185 {
186 // no need to delete the histograms, done because they have been added to the pool
187 }
188 
190 
191 void Coinc2D::DoCanvas(TCanvas *c, Option_t *o)
192 {
193  TString opt = o;
194 
195  if ( opt.Contains("") ) {
196  c->Divide(1,2,0.001,0.001);
197  c->cd(1);
198  fTrackedEnergy->Draw(o);
199  TPad *pad = new TPad("GxG_FOLD", "Pad for Fold distribution",0.7,0.6,1,1);
200  pad->Draw();
201  pad->cd();
202  fFold->Draw();
203  c->cd(2);
204  fTrackedEnergyDopler->Draw(o);
205  }
206 }
207 
209 
210 void Coinc2D::Exec(Option_t */*option*/)
211 {
212 // Frame by Frame action
213  // be sure the frame has been set properly
214  if ( fFrame == 0x0 || !fFrame->IsValid() )
215  { SetLastExecStatus(1); return; }
216 // printf("... is called \n");
217 
218  if ( gDebug > 3 )
219  printf("Coinc2D::::Exec is called \n");
220 
221  // to get the data part of the frame
222  const GammaTrackedInterface *data =
223  GetCstDataPointer<GammaTrackedInterface>(fFrame);
224  VertexInterface *vertex = fVertexBuilder->GetVertex();
225  // vertex->SetDirection(0.,-1.,1.);
226  // vertex->SetPosition(0.,0.,-50.);
227  // vertex->SetBeta(0.02);
228 
229  // now fill histograms
230  UInt_t fold = data->GetNbGamma();
231 
232  fFold->Fill( fold );
233  for (UShort_t i = 0u; i < fold; i++) {
234 
235  const TrackedHit *gamma1 = data->GetGamma(i);
236 
237  fPos3D->Fill(gamma1->GetX(),gamma1->GetY(),gamma1->GetZ());
238 // cout << gamma1->GetHit(1)->GetX() << " " << gamma1->GetHit(1)->GetY() << " " << gamma1->GetHit(1)->GetZ() << endl;
239 
240  Double_t edc1 = DoDopplerCorrection(gamma1, vertex);
241  Double_t e1 = gamma1->GetE();
242 
243  fTrackedEnergy->Fill(e1);
244  if ( vertex->GetBeta() > 0.001 )
245  fTrackedEnergyDopler->Fill(edc1);
246 
247  for (UShort_t j = i+1u; j < fold; j++) {
248 
249  const TrackedHit *gamma2 = data->GetGamma(j);
250 
251  Double_t edc2 = DoDopplerCorrection(gamma2, vertex);
252  Double_t e2 = gamma2->GetE();
253 
254  fGxG->Fill(e1,e2);
255  fGxG->Fill(e2,e1);
256  if ( vertex->GetBeta() > 0.001 ) {
257  fGxG_DC->Fill(edc1,edc2);
258  fGxG_DC->Fill(edc2,edc1);
259  }
260  }
261  }
262 }
263 
267 
268 //Float_t EveTrack::fgMetric = 10.;
269 
270 //EveTrack::EveTrack(const char *name, const char *title, TDirectory *sp_dir, TDirectory *tag_dir):
271 // TrackedWatcher(name, title, sp_dir, tag_dir),
272 // fMaxCycle(100),
273 // fPathFile(Gw::AgataEventDisplay::GetDefaultAgataPath())
274 //{
275 // AgataEventDisplay::Instance()->SetTrackStyle("Cone");
276 
277 // fMeanRadius = MakeTH1<TH1F>("MeanRadius","Mean radius of the tracks (cm)", 500, -1, 34);
278 // fFold = MakeTH1<TH1I>("Fold","Number of Tracks per event", 22, -1, 20) ;
279 // fHitEnergy = MakeTH1<TH1F>("HitEnergy","Energy of hits", 1000,-1,1000.) ;
280 // fTrackEnergy = MakeTH1<TH1F>("TrackEnergy","Energy of tracks", 2000,0.,2000.) ;
281 //}
282 
284 
285 //EveTrack::~EveTrack()
286 //{
287 
288 //}
289 
291 
292 //void EveTrack::SetActive(Bool_t active)
293 //{
294 // TTask::SetActive(active);
295 // AgataEventDisplay::Instance()->Reset();
296 //}
297 
299 
300 //void EveTrack::DoCanvas(TCanvas* c, Option_t* /*o*/)
301 //{
302 // c->Divide(2,2,0.001,0.001);
303 
304 // c->cd(1);
305 // fMeanRadius->Draw();
306 
307 // c->cd(2);
308 // fFold->Draw();
309 
310 // c->cd(3);
311 // fHitEnergy->Draw();
312 
313 // c->cd(4);
314 // fTrackEnergy->Draw();
315 //}
316 
318 
319 //void EveTrack::ShowEve(const Char_t* agataPathFile)
320 //{
321 // AgataEventDisplay::SetDefaultAgataPath(agataPathFile);
322 // AgataEventDisplay::Instance()->BuildDefaultGeometry(false, false);
323 // AgataEventDisplay::Instance()->ShowDisplay();
324 //}
325 
327 
328 //void EveTrack::Exec(Option_t */*option*/)
329 //{
330 // // Frame by Frame action
331 
332 // // be sure the frame has been set properly
333 // // be sure the frame has been set properly
334 // if ( fFrame == 0x0 || !fFrame->IsValid() )
335 // { SetLastExecStatus(1); return; }
336 
337 // if ( gDebug > 3 )
338 // printf("EveTrack::::Exec is called \n");
339 
340 // // to get the data part of the frame
341 // const GammaTrackedInterface *data = GetCstDataPointer<GammaTrackedInterface>(fFrame);
342 
343 // // now fill display container
344 // UInt_t fold = data->GetNbGamma();
345 
346 // if (fold == 0)
347 // return;
348 
349 // fFold->Fill(fold);
350 
351 // AgataEventContainer* agataContainer = AgataEventDisplay::Instance()->GetEventContainer();
352 
353 // if (agataContainer->GetNofEvents() > fMaxCycle) {
354 // TTask::SetActive(false);
355 // }
356 
357 // if ( gDebug > 3 )
358 // std::cout << "TrackHits " << std::endl;
359 
360 // if ( gDebug > 3 )
361 // std::cout << "Fold " << fold << std::endl;
362 
363 // for (UShort_t i = 0u; i < fold; ++i) {
364 
365 // const TrackedHit* gamma = (TrackedHit*)data->GetGamma(i);
366 // TrackHit* hit = agataContainer->NewTrackHit();
367 // fTrackEnergy->Fill(gamma->GetE());
368 
369 // hit->SetX(gamma->GetX()/fgMetric);
370 // hit->SetY(gamma->GetY()/fgMetric);
371 // hit->SetZ(gamma->GetZ()/fgMetric);
372 // hit->SetE(gamma->GetE());
373 
374 // for (Int_t j = 0; j < gamma->GetNbHits(); ++j) {
375 // Hit* gammaHit = gamma->GetHit(j);
376 // StdHit* shit = hit->NewHit();
377 
378 // Float_t radius = TMath::Sqrt(gammaHit->GetX()*gammaHit->GetX() + gammaHit->GetY()*gammaHit->GetY() +
379 // gammaHit->GetZ()*gammaHit->GetZ());
380 // fMeanRadius->Fill(radius/fgMetric);
381 // fHitEnergy->Fill(gammaHit->GetE());
382 // if ( gDebug > 3 ) {
383 // std::cout << "TrackHit Pos " << gammaHit->GetX()/fgMetric << " " << gammaHit->GetY()/fgMetric
384 // << " " << gammaHit->GetZ()/fgMetric << " " << gammaHit->GetE() << std::endl;
385 // }
386 
387 // shit->SetX(gammaHit->GetX()/fgMetric);
388 // shit->SetY(gammaHit->GetY()/fgMetric);
389 // shit->SetZ(gammaHit->GetZ()/fgMetric);
390 // shit->SetE(gammaHit->GetE());
391 // }
392 // }
393 // agataContainer->FillTracks();
394 //}
395 
399 
400 TrackingWatcher::TrackingWatcher(const char *name, const char *title, TDirectory *sp_dir, TDirectory *tag_dir):
401  TrackedWatcher(name,title,sp_dir,tag_dir),
402  CanvasVisu(name,title),
403  fVertexBuilder(0x0),
404  fTrackedEnergy(0x0),
405  fTrackedEnergyDopler(0x0),
406  fGammaMult(0x0),
407  fGxG(0x0),
408  fGxG_DC(0x0),
409  fPos3D(0x0)
410 {
411  // Vertex used for doppler correction, last registered
412  fVertexBuilder = VertexBuilder::theCurrentVertexBuilder();
413 
414  // Now spectra
415  fTrackedEnergy = MakeTH1<TH1F>("TrackedEnergy","Tracked Energy",8192,0,4096);
416 
417  fTrackedEnergyDopler = MakeTH1<TH1F>("TrackedEnergyDC","Tracked Energy with Doppler correction",8192,0,4096);
418  TagOn(fTrackedEnergyDopler);
419 
420  fGammaMult = MakeTH1<TH1F>("GammaMult","Gamma multiplicity",20,0,20);
421  TagOn(fGammaMult);
422 
423  fGxG = MakeTH2<TH2F>("GxG","Gamma Gamma matrix",4096,0,4096,4096,0,4096);
424  fGxG_DC = MakeTH2<TH2F>("GxG_DC","Gamma Gamma matrix with Doppler Correction",4096,0,4096,4096,0,4096);
425 
426  fETheta = MakeTH2<TH2F>("ETheta","Energy vs Theta matrix",4096,0,4096,1800,0,180);
427  fEThetaDC = MakeTH2<TH2F>("ETheta_DC","Energy vs Theta matrix with Doppler Correction",4096,0,4096,1800,0,180);
428 
429  fPos3D = MakeTH3<TH3F>("Pos3D","Positions of Tracked Gamma-rays ;x;y;z",200,-400,400,200,-400,400,200,-400,400);
430 }
431 
433 
435 {
436  // no need to delete the histograms, done because they have been added to the pool
437 }
438 
440 
442 {
443  vertex->SetDirection(0.,-1.,1.);
444  vertex->SetPosition(0.,0.,-50.);
445  vertex->SetBeta(0.02);
446 }
447 
449 
450 void TrackingWatcher::Exec(Option_t */*option*/)
451 {
452  // Frame by Frame action
453  // be sure the frame has been set properly
454  if ( fFrame == 0x0 || !fFrame->IsValid() )
455  { SetLastExecStatus(1); return; }
456 
457  // to get the data part of the frame
458  const GammaTrackedInterface *data =
459  GetCstDataPointer<GammaTrackedInterface>(fFrame);
460  VertexInterface *vertex = fVertexBuilder->GetVertex();
461 
462 
464  // InitVertex(vertex);
465 
466  // now fill histograms
467  UInt_t fold = data->GetNbGamma();
468 
469  fGammaMult->Fill( fold );
470  for (UShort_t i = 0u; i < fold; i++) {
471 
472  const TrackedHit *gamma1 = data->GetGamma(i);
473 
474  double x,y,z;
475  gamma1->GetXYZ(x,y,z);
476 
477  double Theta = GetTheta(gamma1,vertex);
478 
479  fPos3D->Fill(gamma1->GetX(),gamma1->GetY(),gamma1->GetZ());
480 
481  Double_t edc1 = DoDopplerCorrection(gamma1, vertex);
482  Double_t e1 = gamma1->GetE();
483 
484  fTrackedEnergy->Fill(e1);
485  fETheta->Fill(e1,Theta);
486  fEThetaDC->Fill(edc1,Theta);
487 
488  if ( vertex->GetBeta() > 0.001 )
489  fTrackedEnergyDopler->Fill(edc1);
490 
491  for (UShort_t j = i+1u; j < fold; j++) {
492 
493  const TrackedHit *gamma2 = data->GetGamma(j);
494 
495  Double_t edc2 = DoDopplerCorrection(gamma2, vertex);
496  Double_t e2 = gamma2->GetE();
497 
498  fGxG->Fill(e1,e2);
499  fGxG->Fill(e2,e1);
500  if ( vertex->GetBeta() > 0.001 ) {
501  fGxG_DC->Fill(edc1,edc2);
502  fGxG_DC->Fill(edc2,edc1);
503  }
504  }
505  }
506 }
507 
509 
511 {
512  Double_t px, py, pz, dx, dy, dz;
513  vertex->GetPosition (px, py, pz);
514  vertex->GetDirection(dx, dy, dz);
515 
516  Double_t x = hit->GetX() - px;
517  Double_t y = hit->GetY() - py;
518  Double_t z = hit->GetZ() + fZOffset - pz;
519 
520  Double_t cosTheta=0;
521 
522  Double_t dd = x*x + y*y + z*z;
523 
524  if( !(dd == 0.0) ) cosTheta = (x*dx + y*dy + z*dz)/::sqrt(dd);
525 
526  return TMath::ACos(cosTheta)*TMath::RadToDeg();
527 }
528 
530 
532 {
533  TString Name = "Tracked_Energy";
534  TString Title = "tracked energy without and with doppler correction";
535 
536  TCanvas *c = CanvasVisu::NewCanvas(Name,Title);
537 
538  TH1 *h[2] = {fTrackedEnergy,fTrackedEnergyDopler};
539 
540  if(!OnSamePad)
541  {
542  c->Divide(1,2,0.001,0.001);
543 
544  double topmarg = 0.0657599;
545  double rightmarg = 0.00726832;
546  double leftmarg = 0.0490612;
547  double bottommarg = 0.114119;
548 
549  for(int i=0 ; i<2 ; i++)
550  {
551  c->cd(i+1);
552 
553  gPad->SetRightMargin(rightmarg);
554  gPad->SetTopMargin(topmarg);
555  gPad->SetLeftMargin(leftmarg);
556  gPad->SetBottomMargin(bottommarg);
557 
558  TH1 *htemp = h[i];
559 
560  htemp->SetLineColor(GetColor(0));
561  htemp->Draw();
562 
563  htemp->SetStats(false);
564  htemp->GetXaxis()->SetTitle("Energy (keV)");
565  htemp->GetXaxis()->SetTitleSize(0.06);
566  htemp->GetXaxis()->SetLabelSize(0.06);
567  htemp->GetXaxis()->SetTitleOffset(0.91);
568 
569  htemp->GetYaxis()->SetTitle("Counts");
570  htemp->GetYaxis()->SetTitleSize(0.06);
571  htemp->GetYaxis()->SetLabelSize(0.06);
572  htemp->GetYaxis()->SetTitleOffset(0.45);
573 
574  gPad->Modified();
575  gPad->Update();
576 
577  TPaveStats *st = (TPaveStats*)htemp->FindObject("stats");
578  if(st)
579  {
580  st->SetX2NDC(1-rightmarg);
581  st->SetY2NDC(1-topmarg);
582  }
583 
584  TPaveText *txt = (TPaveText*)gPad->FindObject("title");
585  if(txt)
586  {
587  double width = txt->GetY2NDC()-txt->GetY1NDC();
588  txt->SetY2NDC(1-topmarg);
589  txt->SetY1NDC(txt->GetY2NDC()-width);
590  }
591 
592  if(i==0)
593  {
594  TPad *pad = new TPad("GammaMult", "Gamma multiplicity",0.7,0.6,1-rightmarg-0.0025,1-topmarg-0.005);
595  pad->SetRightMargin(0.00486268);
596  pad->SetTopMargin(0.00912681);
597  pad->SetLeftMargin(0.0736404);
598  pad->SetBottomMargin(0.195208);
599  pad->Draw();
600  pad->cd();
601  pad->SetLogy();
602  fGammaMult->SetLineColor(GetColor(0));
603  fGammaMult->SetStats(false);
604  fGammaMult->GetXaxis()->SetTitle("Gamma multiplicity");
605  fGammaMult->GetXaxis()->SetTitleSize(0.08);
606  fGammaMult->GetXaxis()->SetLabelSize(0.08);
607 
608  fGammaMult->GetYaxis()->SetTitle("Counts");
609  fGammaMult->GetYaxis()->SetTitleSize(0.08);
610  fGammaMult->GetYaxis()->SetLabelSize(0.08);
611  fGammaMult->GetYaxis()->SetTitleOffset(0.5);
612  fGammaMult->Draw();
613  }
614  }
615  }
616  else
617  {
618  double topmarg = 0.00240964;
619  double rightmarg = 0.0169287;
620  double leftmarg = 0.0719468;
621  double bottommarg = 0.0843373;
622 
623  c->cd();
624 
625  gPad->SetRightMargin(rightmarg);
626  gPad->SetTopMargin(topmarg);
627  gPad->SetLeftMargin(leftmarg);
628  gPad->SetBottomMargin(bottommarg);
629 
630  h[0]->SetBit(TH1::kNoTitle);
631  h[0]->SetLineColor(GetColor(0));
632  h[0]->SetStats(false);
633  h[0]->Draw();
634 
635  h[0]->GetXaxis()->SetTitle("Energy (keV)");
636  h[0]->GetXaxis()->SetTitleSize(0.04);
637  h[0]->GetXaxis()->SetLabelSize(0.04);
638  h[0]->GetXaxis()->SetTitleOffset(1.);
639 
640  h[0]->GetYaxis()->SetTitle("Counts");
641  h[0]->GetYaxis()->SetTitleSize(0.04);
642  h[0]->GetYaxis()->SetLabelSize(0.04);
643  h[0]->GetYaxis()->SetTitleOffset(1.);
644 
645  h[1]->SetLineColor(GetColor(1));
646  h[1]->SetStats(false);
647  h[1]->Draw("same");
648 
649  gPad->Modified();
650  gPad->Update();
651 
652  TLegend *leg = new TLegend(0.756348,0.870486,1-rightmarg,1-topmarg,"Tracked Energy","NDC");
653  leg->AddEntry(h[0],"Without doppler correction","l");
654  leg->AddEntry(h[1],"With doppler correction","l");
655  leg->Draw();
656  }
657 }
658 
660 
661 void TrackingWatcher::ShowGammaGamma(const char *NoDC_DC_Both)
662 {
663  TString Name = "Tracking_GammaGamma";
664  TString Title = "Gamma Gamma matrix of the tracked gamma";
665 
666  TCanvas *c = CanvasVisu::NewCanvas(Name,Title);
667 
668  if(((TString)NoDC_DC_Both) == "Both")
669  {
670  c->Divide(1,2,0.001,0.001);
671 
672  double topmarg = 0.0318624;
673  double rightmarg = 0.00726832;
674  double leftmarg = 0.0529148;
675  double bottommarg = 0.114119;
676 
677  TH2 *h[2] = {fGxG,fGxG_DC};
678 
679  for(int i=0 ; i<2 ; i++)
680  {
681  c->cd(i+1);
682 
683  gPad->SetRightMargin(rightmarg);
684  gPad->SetTopMargin(topmarg);
685  gPad->SetLeftMargin(leftmarg);
686  gPad->SetBottomMargin(bottommarg);
687 
688  TH1 *htemp = h[i];
689 
690  htemp->SetLineColor(GetColor(0));
691  htemp->Draw("col");
692  gPad->SetLogz();
693 
694  htemp->SetBit(TH2::kNoTitle);
695  htemp->GetXaxis()->SetTitle("Energy #gamma_{1} (keV)");
696  htemp->GetXaxis()->SetTitleSize(0.06);
697  htemp->GetXaxis()->SetLabelSize(0.06);
698  htemp->GetXaxis()->SetTitleOffset(0.91);
699 
700  htemp->GetYaxis()->SetTitle("Energy #gamma_{2} (keV)");
701  htemp->GetYaxis()->SetTitleSize(0.06);
702  htemp->GetYaxis()->SetLabelSize(0.06);
703  htemp->GetYaxis()->SetTitleOffset(0.45);
704 
705  gPad->Modified();
706  gPad->Update();
707 
708  TPaveStats *st = (TPaveStats*)htemp->FindObject("stats");
709  if(st)
710  {
711  st->SetX2NDC(1-rightmarg);
712  st->SetY2NDC(1-topmarg);
713  }
714 
715  TPaveText *txt = (TPaveText*)gPad->FindObject("title");
716  if(txt)
717  {
718  double width = txt->GetY2NDC()-txt->GetY1NDC();
719  txt->SetY2NDC(1-topmarg);
720  txt->SetY1NDC(txt->GetY2NDC()-width);
721  }
722  }
723  }
724  else
725  {
726  TH2 *htemp = 0x0;
727  if(((TString)NoDC_DC_Both) == "NoDC") htemp = fGxG;
728  else if(((TString)NoDC_DC_Both) == "DC") htemp = fGxG_DC;
729 
730  double topmarg = 0.0156627;
731  double rightmarg = 0.0169287;
732  double leftmarg = 0.0719468;
733  double bottommarg = 0.0843373;
734 
735  gPad->SetRightMargin(rightmarg);
736  gPad->SetTopMargin(topmarg);
737  gPad->SetLeftMargin(leftmarg);
738  gPad->SetBottomMargin(bottommarg);
739 
740  htemp->SetBit(TH2::kNoTitle);
741  htemp->GetXaxis()->SetTitle("Energy #gamma_{1} (keV)");
742  htemp->GetYaxis()->SetTitle("Energy #gamma_{2} (keV)");
743  htemp->Draw("col");
744  gPad->SetLogz();
745 
746  gPad->Modified();
747  gPad->Update();
748 
749  TPaveStats *st = (TPaveStats*)htemp->FindObject("stats");
750  if(st)
751  {
752  st->SetX2NDC(1-rightmarg);
753  st->SetY2NDC(1-topmarg);
754  }
755  }
756 }
757 
759 
760 void TrackingWatcher::ShowETheta(const char *NoDC_DC_Both)
761 {
762  TString Name = "ETheta";
763  TString Title = "Energy Theta matrix of the tracked gamma";
764 
765  TCanvas *c = CanvasVisu::NewCanvas(Name,Title);
766 
767  TH2 *h[2] = {fETheta,fEThetaDC};
768 
769  if(((TString)NoDC_DC_Both) == "Both")
770  {
771  c->Divide(1,2,0.001,0.001);
772 
773  double topmarg = 0.0318624;
774  double rightmarg = 0.00726832;
775  double leftmarg = 0.0529148;
776  double bottommarg = 0.114119;
777 
778  for(int i=0 ; i<2 ; i++)
779  {
780  c->cd(i+1);
781 
782  gPad->SetRightMargin(rightmarg);
783  gPad->SetTopMargin(topmarg);
784  gPad->SetLeftMargin(leftmarg);
785  gPad->SetBottomMargin(bottommarg);
786 
787  TH1 *htemp = h[i];
788 
789  htemp->SetLineColor(GetColor(0));
790  htemp->Draw("col");
791  gPad->SetLogz();
792 
793  htemp->SetBit(TH2::kNoTitle);
794  htemp->GetXaxis()->SetTitle("Energy (keV)");
795  htemp->GetXaxis()->SetTitleSize(0.06);
796  htemp->GetXaxis()->SetLabelSize(0.06);
797  htemp->GetXaxis()->SetTitleOffset(0.91);
798 
799  htemp->GetYaxis()->SetTitle("#theta (deg)");
800  htemp->GetYaxis()->SetTitleSize(0.06);
801  htemp->GetYaxis()->SetLabelSize(0.06);
802  htemp->GetYaxis()->SetTitleOffset(0.45);
803 
804  gPad->Modified();
805  gPad->Update();
806 
807  TPaveStats *st = (TPaveStats*)htemp->FindObject("stats");
808  if(st)
809  {
810  st->SetX2NDC(1-rightmarg);
811  st->SetY2NDC(1-topmarg);
812  }
813 
814  TPaveText *txt = (TPaveText*)gPad->FindObject("title");
815  if(txt)
816  {
817  double width = txt->GetY2NDC()-txt->GetY1NDC();
818  txt->SetY2NDC(1-topmarg);
819  txt->SetY1NDC(txt->GetY2NDC()-width);
820  }
821  }
822  }
823  else
824  {
825  TH2 *htemp = 0x0;
826  if(((TString)NoDC_DC_Both) == "NoDC") htemp = h[0];
827  else if(((TString)NoDC_DC_Both) == "DC") htemp = h[1];
828 
829  double topmarg = 0.0156627;
830  double rightmarg = 0.0169287;
831  double leftmarg = 0.0719468;
832  double bottommarg = 0.0843373;
833 
834  gPad->SetRightMargin(rightmarg);
835  gPad->SetTopMargin(topmarg);
836  gPad->SetLeftMargin(leftmarg);
837  gPad->SetBottomMargin(bottommarg);
838 
839  htemp->SetBit(TH2::kNoTitle);
840  htemp->GetXaxis()->SetTitle("Energy #gamma_{1} (keV)");
841  htemp->GetYaxis()->SetTitle("#theta (deg))");
842  htemp->Draw("col");
843  gPad->SetLogz();
844 
845  gPad->Modified();
846  gPad->Update();
847 
848  TPaveStats *st = (TPaveStats*)htemp->FindObject("stats");
849  if(st)
850  {
851  st->SetX2NDC(1-rightmarg);
852  st->SetY2NDC(1-topmarg);
853  }
854  }
855 }
856 
858 
859 void TrackingWatcher::ShowCommonSumTrackedDC(bool OnSamePad, const char *DC_NoDC)
860 {
861  TH1 *CoreCommon = (TH1*)gROOT->FindObjectAny("CoreCommon");
862  TH1 *SumSpectra = (TH1*)gROOT->FindObjectAny("SumSpectra");
863 
864  if(!CoreCommon)
865  {
866  cout<<"CoreCommon spectrum not found"<<endl;
867  return;
868  }
869  if(!SumSpectra)
870  {
871  cout<<"SumSpectra spectrum not found"<<endl;
872  return;
873  }
874 
875  TH1 *trackedhist = 0x0;
876  if(((TString)DC_NoDC) == "DC") trackedhist = fTrackedEnergyDopler;
877  else if(((TString)DC_NoDC) == "NoDC") trackedhist = fTrackedEnergy;
878  else
879  {
880  cout<<"Wrong option, DC_NoDC should be DC or NoDC"<<endl;
881  return;
882  }
883 
884  TString Name = "CommonSumTrackedDC";
885  TString Title = "comparison of Core common, sum spectra and tracked spectrum dc";
886 
887  TCanvas *c = CanvasVisu::NewCanvas(Name,Title);
888 
889  TH1 *h[3] = {CoreCommon,SumSpectra,trackedhist};
890 
891  if(!OnSamePad)
892  {
893  c->Divide(1,3,0.001,0.001);
894 
895  double topmarg = 0.0657599;
896  double rightmarg = 0.00142122;
897  double leftmarg = 0.0444335;
898  double bottommarg = 0.127278;
899 
900  for(int i=0 ; i<3 ; i++)
901  {
902  c->cd(i+1);
903 
904  gPad->SetRightMargin(rightmarg);
905  gPad->SetTopMargin(topmarg);
906  gPad->SetLeftMargin(leftmarg);
907  gPad->SetBottomMargin(bottommarg);
908 
909  TH1 *htemp = h[i];
910 
911  htemp->SetLineColor(GetColor(0));
912  htemp->Draw();
913 
914  htemp->SetStats(false);
915  htemp->GetXaxis()->SetTitle("Energy (keV)");
916  htemp->GetXaxis()->SetTitleSize(0.07);
917  htemp->GetXaxis()->SetLabelSize(0.07);
918  htemp->GetXaxis()->SetTitleOffset(0.87);
919 
920  htemp->GetYaxis()->SetTitle("Counts");
921  htemp->GetYaxis()->SetTitleSize(0.07);
922  htemp->GetYaxis()->SetLabelSize(0.07);
923  htemp->GetYaxis()->SetTitleOffset(0.35);
924 
925  gPad->Modified();
926  gPad->Update();
927 
928  TPaveStats *st = (TPaveStats*)htemp->FindObject("stats");
929  if(st)
930  {
931  st->SetX2NDC(1-rightmarg);
932  st->SetY2NDC(1-topmarg);
933  }
934 
935  TPaveText *txt = (TPaveText*)gPad->FindObject("title");
936  if(txt)
937  {
938  double width = txt->GetY2NDC()-txt->GetY1NDC();
939  txt->SetY2NDC(1-topmarg);
940  txt->SetY1NDC(txt->GetY2NDC()-width);
941  }
942  }
943  }
944  else
945  {
946  double topmarg = 0.00240964;
947  double rightmarg = 0.0169287;
948  double leftmarg = 0.0719468;
949  double bottommarg = 0.0843373;
950 
951  c->cd();
952 
953  gPad->SetRightMargin(rightmarg);
954  gPad->SetTopMargin(topmarg);
955  gPad->SetLeftMargin(leftmarg);
956  gPad->SetBottomMargin(bottommarg);
957 
958  h[0]->SetBit(TH1::kNoTitle);
959  h[0]->SetLineColor(GetColor(0));
960  h[0]->SetStats(false);
961  h[0]->Draw();
962 
963  h[0]->GetXaxis()->SetTitle("Energy (keV)");
964  h[0]->GetXaxis()->SetTitleSize(0.04);
965  h[0]->GetXaxis()->SetLabelSize(0.04);
966  h[0]->GetXaxis()->SetTitleOffset(1.);
967 
968  h[0]->GetYaxis()->SetTitle("Counts");
969  h[0]->GetYaxis()->SetTitleSize(0.04);
970  h[0]->GetYaxis()->SetLabelSize(0.04);
971  h[0]->GetYaxis()->SetTitleOffset(1.);
972 
973  h[1]->SetLineColor(GetColor(1));
974  h[1]->SetStats(false);
975  h[1]->Draw("same");
976 
977  h[2]->SetLineColor(GetColor(2));
978  h[2]->SetStats(false);
979  h[2]->Draw("same");
980 
981  gPad->Modified();
982  gPad->Update();
983 
984  TLegend *leg = new TLegend(0.706167,0.870486,1-rightmarg,1-topmarg,"Core common -vs- Sum spectra (Calorimeter) -vs- tracked energy","NDC");
985  leg->AddEntry(h[0],"Core common","l");
986  leg->AddEntry(h[1],"Sum spectra (Calorimeter)","l");
987 
988  if(((TString)DC_NoDC) == "DC") leg->AddEntry(h[2],"Tracked, doppler corrected","l");
989  else leg->AddEntry(h[2],"Tracked, no doppler corrected","l");
990 
991  leg->Draw();
992  }
993 }
994 
996 
998 {
999  TString Name = "Tracked_3DHits";
1000  TString Title = "3D view of the first hit of the tracked gamma";
1001 
1002  CanvasVisu::NewCanvas(Name,Title);
1003 
1004  fPos3D->Draw();
1005 
1006 }
1007 
1009 ClassImp(Coinc2D);
1010 //ClassImp(EveTrack);
virtual Bool_t IsValid() const
true if it is a valid pointer
Definition: Frame.h:616
virtual void GetDirection(Double_t &, Double_t &, Double_t &, Double_t=0.0) const =0
get the direction of the source (last argument is used in case the position depends on time) ...
virtual VertexInterface * GetVertex()
Get the vertex data interface.
virtual Double_t DoDopplerCorrection(const TrackedHit *hit, VertexInterface *vertex)
Do Doppler taking into account additionnal offset of the agat position (see SetAgataOffset) ...
printf("******************************************************************** \n")
virtual Double_t GetY() const
Definition: Hits.h:57
virtual void Exec(Option_t *option="")
watch the current frame
virtual ~TrackingWatcher()
******************************************************************************************/// ...
void ShowCommonSumTrackedDC(bool OnSamePad=false, const char *DC_NoDC="DC")
******************************************************************************************/// ...
virtual void GetPosition(Double_t &, Double_t &, Double_t &, Double_t=0.0) const =0
get the position of the source (last argument is used in case the position depends on time) ...
virtual Double_t GetE() const
Definition: Hits.h:67
virtual void DoCanvas(TCanvas *c, Option_t *)
To be overwritten by real implementation if a canvas is produced.
The tracking algorithm produces a stack of TrackedHits.
Definition: TrackedFrame.h:41
virtual void SetPosition(Double_t, Double_t, Double_t, Double_t=0.0)=0
Set the position of the source (last argument is used in case the position depends on time) ...
TH1F * fTrackedEnergyDopler
virtual Bool_t SetTrigger(ADF::DFTrigger *=0x0)
Set the trigger attached to this watcher.
virtual TrackedHit * GetGamma(UShort_t) const =0
To get the current number of gammas in the stack.
TCanvas * NewCanvas(TString cname, TString ctitle)
Definition: CanvasVisu.C:102
virtual Double_t GetBeta(Double_t=0.0) const =0
get recoil velocity
void Show3DHits()
******************************************************************************************/// ...
double GetTheta(const TrackedHit *hit, VertexInterface *vertex)
******************************************************************************************/// ...
void SetLastExecStatus(Short_t s=0)
reset last status. 0 means no error, 0 < means error, > 0 means ok with conditions ...
Definition: Watchers.h:294
virtual Double_t GetT() const
Definition: Hits.h:72
virtual ~Coinc2D()
******************************************************************************************/// ...
void ShowTrackedSpectra(bool OnSamePad=false)
******************************************************************************************/// ...
virtual void DopplerCorrection(Hit *)=0
header file for AgataKeyFactory.cpp
void ShowETheta(const char *NoDC_DC_Both="DC")
******************************************************************************************/// ...
It is a hit associated to a list of Hits.
Definition: Hits.h:256
Color_t GetColor(int i)
Definition: CanvasVisu.h:49
virtual void SetAgataOffset(Double_t offset=0.0)
set the agata offset in case it is not @ the optimal position
virtual Double_t GetX() const
Definition: Hits.h:55
TH1F * fTrackedEnergy
void TagOn(TObject *)
Add this histogram to the list of tagged histograms.
Definition: Watchers.cpp:901
virtual Double_t GetGamma(Double_t=0.0) const
Definition: MetaFrame.h:77
TrackingWatcher(const char *name, const char *title, TDirectory *sp_dir=0x0, TDirectory *tag_dir=0x0)
*/
Base class for a trigger on a data flow.
Definition: Trigger.h:155
ADF::LogMessage & endl(ADF::LogMessage &log)
void InitVertex(VertexInterface *vertex)
******************************************************************************************/// ...
virtual void SetDirection(Double_t, Double_t, Double_t, Double_t=0.0)=0
set the direction of the source (last argument is used in case the position depends on time) ...
virtual void GetXYZ(Double_t &x, Double_t &y, Double_t &z) const
Definition: Hits.h:64
virtual Double_t GetZ() const
Definition: Hits.h:59
virtual UShort_t GetNbGamma() const =0
To get the current number of gammas in the stack.
virtual void SetBeta(Double_t, Double_t=0.0)=0
set recoil velocity
void ShowGammaGamma(const char *NoDC_DC_Both="DC")
******************************************************************************************/// ...
Coinc2D(const char *name, const char *title, TDirectory *sp_dir=0x0, TDirectory *tag_dir=0x0)
*/
LogMessage & endl(LogMessage &)
Base Watcher working for any kind of Trackek Frame (Frame interface) and.
ClassImp(TrackedWatcher)
virtual void Exec(Option_t *option="")
watch the current frame
static VertexBuilder * theCurrentVertexBuilder()
to get the current VertexWatcher (for other watchers) i.e. the last one registered.