R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandGainMatching.h
Go to the documentation of this file.
1/******************************************************************************
2 * Copyright (C) 2019 GSI Helmholtzzentrum für Schwerionenforschung GmbH *
3 * Copyright (C) 2019-2025 Members of R3B Collaboration *
4 * *
5 * This software is distributed under the terms of the *
6 * GNU General Public Licence (GPL) version 3, *
7 * copied verbatim in the file "LICENSE". *
8 * *
9 * In applying this license GSI does not waive the privileges and immunities *
10 * granted to it by virtue of its status as an Intergovernmental Organization *
11 * or submit itself to any jurisdiction. *
12 ******************************************************************************/
13
14#ifndef R3BNEULANDGAINMATCHING_H
15#define R3BNEULANDGAINMATCHING_H
16
17#include "FairTask.h"
18#include "R3BChannelAccessEPICS.h"
19#include "R3BEventHeader.h"
20#include "TH1.h"
21#include <TClonesArray.h>
22
23class R3BNeulandGainMatching : public FairTask
24{
25
26 public:
32
39 R3BNeulandGainMatching(const char* name, Int_t iVerbose = 1);
40
46
53 virtual InitStatus Init();
54
60 virtual void Exec(Option_t* option);
61
67 virtual void FinishEvent();
68
73 virtual void FinishTask();
74
79 inline void SetUpdateRate(Int_t rate) { fUpdateRate = rate; }
80
85 inline void SetTrigger(Int_t trigger) { fTrigger = trigger; }
86
91 inline void SetNofModules(Int_t firstPlane, Int_t nPlanes, Int_t nPaddles)
92 {
93 fFirstPlane = firstPlane;
94 fNofPlanes = nPlanes;
95 fNofBarsPerPlane = nPaddles;
96 fNofPMTs = nPlanes * nPaddles * 2;
97 }
98
102 inline void SetNeededStat(Int_t nevents) { fNEventsNeeded = nevents; }
103
104 private:
105 UInt_t fNofPlanes;
106 UInt_t fNofBarsPerPlane;
107 UInt_t fNofPMTs;
108
109 UInt_t fFirstPlane;
110
111 Int_t fUpdateRate;
112 Int_t fTrigger;
113 Int_t fNEvents;
114
115 Int_t fNEventsNeeded;
116
117 Bool_t finished;
118
119 R3BChannelAccessMasterEPICS epics;
120
121 struct
122 {
123 R3BChannelAccessGroup* group;
124 R3BChannelAccess* vmon;
125 R3BChannelAccess* vtarget;
126 } ca[60][50][2];
127
128 TH1F* hCosmicPeak[60][50][2];
129
130 Int_t iteration[60][50][2];
131 Double_t esum[60][50][2];
132 Double_t ealt[60][50][2];
133
134 Double_t hv[60][50][2];
135
136 TClonesArray* fPmt;
137 R3BEventHeader* header;
138
139 Int_t peakmethod;
140 Double_t getcosmicpeak(TH1*);
141 Double_t searchcosmicpeak(TH1*, Double_t);
142 Double_t searchcosmicpeaknb(TH1*, Double_t);
143
144 public:
145 ClassDef(R3BNeulandGainMatching, 1)
146};
147
148#endif
virtual void Exec(Option_t *option)
Method for event loop implementation.
void SetNeededStat(Int_t nevents)
Method for setting needed statistics for data analyis.
virtual InitStatus Init()
Method for task initialization.
void SetUpdateRate(Int_t rate)
Method for setting the update rate.
virtual ~R3BNeulandGainMatching()
Destructor.
void SetTrigger(Int_t trigger)
Method for selecting events with certain trigger value.
R3BNeulandGainMatching()
Default constructor.
void SetNofModules(Int_t firstPlane, Int_t nPlanes, Int_t nPaddles)
Method for setting number of modules in NeuLAND setup.
virtual void FinishEvent()
A method for finish of processing of an event.
virtual void FinishTask()
Method for finish of the task execution.
R3BChannelAccessGroup * group
constexpr Int_t nPlanes