CbmRoot
Loading...
Searching...
No Matches
CbmRichConverter.h
Go to the documentation of this file.
1/* Copyright (C) 2012-2019 UGiessen/JINR-LIT, Giessen/Dubna
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Semen Lebedev [committer] */
4
13#ifndef CBM_RICH_CONVERTER
14#define CBM_RICH_CONVERTER
15
16#include "CbmRichHit.h"
17#include "CbmRichRing.h"
18#include "CbmRichRingLight.h"
19#include "FairRootManager.h"
20#include "TClonesArray.h"
21
22#include <iostream>
23#include <vector>
24
25using std::cout;
26using std::endl;
27using std::vector;
28
39 public:
45 static void CopyHitsToRingLight(const CbmRichRing* ring1, CbmRichRingLight* ring2)
46 {
47 if (NULL == fRichHits) {
48 Init();
49 }
50 if (NULL == fRichHits) return;
51 int nofHits = ring1->GetNofHits();
52 for (int i = 0; i < nofHits; i++) {
53 Int_t hitInd = ring1->GetHit(i);
54 CbmRichHit* hit = (CbmRichHit*) fRichHits->At(hitInd);
55 if (NULL == hit) continue;
56 CbmRichHitLight hl(hit->GetX(), hit->GetY(), hitInd);
57 ring2->AddHit(hl);
58 }
59 }
60
67 static void CopyHitsToRingLight(const vector<double>& hitX, const vector<double>& hitY, CbmRichRingLight* ring)
68 {
69 int nofHits = hitX.size();
70 for (int i = 0; i < nofHits; i++) {
71 CbmRichHitLight hl(hitX[i], hitY[i]);
72 ring->AddHit(hl);
73 }
74 }
75
81 static void CopyParamsToRing(const CbmRichRingLight* ring1, CbmRichRing* ring2)
82 {
83 ring2->SetCenterX(ring1->GetCenterX());
84 ring2->SetCenterY(ring1->GetCenterY());
85 ring2->SetChi2(ring1->GetChi2());
86 ring2->SetAaxis(ring1->GetAaxis());
87 ring2->SetBaxis(ring1->GetBaxis());
88 ring2->SetRadius(ring1->GetRadius());
89 ring2->SetPhi(ring1->GetPhi());
90 }
91
92 static TClonesArray* fRichHits;
93
94 public:
98 static void Init()
99 {
100 FairRootManager* ioman = FairRootManager::Instance();
101 if (NULL == ioman) {
102 cout << "-E- CbmRichConverter::Init, RootManager not instantised!" << endl;
103 return;
104 }
105 fRichHits = (TClonesArray*) ioman->GetObject("RichHit");
106 if (NULL == fRichHits) {
107 cout << "-W- CbmRichConverter::Init, No RichHit array" << endl;
108 }
109 }
110};
111
112
113#endif
double GetY() const
Definition CbmPixelHit.h:74
double GetX() const
Definition CbmPixelHit.h:73
Convert internal data classes to cbmroot common data classes.
static void CopyParamsToRing(const CbmRichRingLight *ring1, CbmRichRing *ring2)
Copy parameters from CbmRichRingLight to CbmRichRing.
static void Init()
Initialize array of RICH hits.
static void CopyHitsToRingLight(const vector< double > &hitX, const vector< double > &hitY, CbmRichRingLight *ring)
Copy hits coordinates from vectors to CbmRichRingLight.
static TClonesArray * fRichHits
static void CopyHitsToRingLight(const CbmRichRing *ring1, CbmRichRingLight *ring2)
Copy hits from CbmRichRing to CbmRichRingLight.
float GetCenterX() const
float GetAaxis() const
float GetRadius() const
float GetCenterY() const
float GetPhi() const
float GetBaxis() const
float GetChi2() const
void AddHit(CbmRichHitLight hit)
Add new hit to the ring.
void SetChi2(double chi2)
Definition CbmRichRing.h:61
uint32_t GetHit(int32_t i) const
Definition CbmRichRing.h:39
void SetAaxis(double a)
Definition CbmRichRing.h:55
int32_t GetNofHits() const
Definition CbmRichRing.h:37
void SetPhi(double phi)
Definition CbmRichRing.h:60
void SetCenterY(float y)
Definition CbmRichRing.h:53
void SetCenterX(float x)
Definition CbmRichRing.h:52
void SetRadius(float r)
Definition CbmRichRing.h:54
void SetBaxis(double b)
Definition CbmRichRing.h:56