Jpp 20.0.0-72-g597b30bc9
the software that should make you happy
Loading...
Searching...
No Matches
JOMGsim.cc File Reference

Auxiliary program to determine PMT parameters from OMGSim data. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <limits>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "TRandom3.h"
#include "km3net-dataformat/offline/Head.hh"
#include "km3net-dataformat/offline/Evt.hh"
#include "km3net-dataformat/online/JDAQ.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JPMTParameters.hh"
#include "JDetector/JPMTAnalogueSignalProcessor.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JMonteCarloFileSupportkit.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JMeta.hh"
#include "JTools/JCombinatorics.hh"
#include "JTools/JRouter.hh"
#include "JAAnet/JHead.hh"
#include "JCalibrate/JCalibrateK40.hh"
#include "JROOT/JRootToolkit.hh"
#include "JROOT/JRootFileWriter.hh"
#include "Jeep/JPrint.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Detailed Description

Auxiliary program to determine PMT parameters from OMGSim data.

Author
mdejong

Definition in file JOMGsim.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 49 of file JOMGsim.cc.

50{
51 using namespace std;
52 using namespace JPP;
53 using namespace KM3NETDAQ;
54
56 JLimit_t& numberOfEvents = inputFile.getLimit();
57 string outputFile;
58 string detectorFile;
59 Double_t TMax_ns;
60 double R;
61 JPMTParameters parameters;
62 double QE = 1.0;
63 ULong_t seed;
64 int debug;
65
66 parameters.TTS_ns = 0.0;
67 parameters.gainSpread = 0.0;
68 parameters.PunderAmplified = 0.0;
69 parameters.slewing = false;
70
71 try {
72
73 JProperties properties = parameters.getProperties();
74
75 properties["QE"] = QE;
76
77 JParser<> zap("Auxiliary program to determine PMT parameters from OMGsim data.");
78
79 zap['f'] = make_field(inputFile, "input file.");
80 zap['o'] = make_field(outputFile, "output file.") = "omgsim.root";
81 zap['n'] = make_field(numberOfEvents) = JLimit::max();
82 zap['a'] = make_field(detectorFile, "detector file.");
83 zap['T'] = make_field(TMax_ns, "time window [ns].") = 20.0;
84 zap['R'] = make_field(R, "radioactivity [Bq]");
85 zap['P'] = make_field(properties, "PMT parameters") = JPARSER::initialised();
86 zap['S'] = make_field(seed, "seed") = 0;
87 zap['d'] = make_field(debug, "debug flag.") = 1;
88
89 zap(argc, argv);
90 }
91 catch(const exception &error) {
92 FATAL(error.what() << endl);
93 }
94
95
96 gRandom->SetSeed(seed);
97
98
100
101 try {
102 load(detectorFile, detector);
103 }
104 catch(const JException& error) {
105 FATAL(error);
106 }
107
108 if (detector.size() != 1u) {
109 FATAL("Wrong number of modules " << detector.size() << endl);
110 }
111
112
113 const JModule& module = detector[0];
114
115 JCombinatorics combinatorics(module.size());
116
117 combinatorics.sort(JPairwiseComparator(module));
118
119 const Double_t WS = 0.25;
120
121 const Int_t nx = combinatorics.getNumberOfPairs();
122 const Double_t xmin = -0.5;
123 const Double_t xmax = nx - 0.5;
124
125 const Double_t ymin = -floor(TMax_ns) + 0.5*WS;
126 const Double_t ymax = +floor(TMax_ns) - 0.5*WS;
127 const Int_t ny = (Int_t) ((ymax - ymin) / WS);
128
129 TH2D h2r(MAKE_CSTRING(module.getID() << _2R), NULL, nx, xmin, xmax, ny, ymin, ymax);
130
131 h2r.Sumw2();
132
133
134 JRouter<int> router;
135
136 for (size_t i = 0; i != module.size(); ++i) {
137 router.put(module[i].getID(), i);
138 }
139
140 double numberOfPrimaries = 0.0;
141
142 Head header;
143
144 try {
145
146 STATUS("Extracting header data... " << flush);
147
148 header = getHeader(inputFile);
149
150 const JHead buffer(header);
151
152 if (buffer.is_valid(&JHead::norma))
153 numberOfPrimaries = buffer.norma.numberOfPrimaries;
154 else
155 THROW(JException, "Missing header information.");
156
157 STATUS("OK" << endl);
158
159 } catch(const exception& error) {
160 FATAL(error.what() << endl);
161 }
162
163 if (numberOfPrimaries == 0.0) {
164 FATAL("Number of primaries " << numberOfPrimaries << endl);
165 }
166
167
168 parameters.QE *= QE;
169
170 const JPMTAnalogueSignalProcessor cpu(parameters);
171
173 JPMTData<JPMTPulse> output;
174
175 struct hit_type {
176 int pmt;
177 double t1;
178 };
179
180 vector<hit_type> buffer;
181
182
183 TFile out(outputFile.c_str(), "recreate");
184
185 putObject(out, JMeta(argc, argv));
186
187 while (inputFile.hasNext()) {
188
189 if (inputFile.getCounter()%10000== 0) {
190 STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
191 }
192
193 const Evt* evt = inputFile.next();
194
195 if (evt->mc_hits.size() >= 2u) {
196
197 buffer.clear();
198
199 for (const auto& hit : evt->mc_hits) {
200
201 if (router.has(hit.pmt_id)) {
202
203 const int pmt = router.get(hit.pmt_id);
204
205 input .clear();
206 output.clear();
207
208 input.push_back(JPMTSignal(hit.t, hit.a));
209
210 cpu(module[pmt].getCalibration(), input, output);
211
212 for (const auto& hit : output) {
213 buffer.push_back({pmt, hit.t_ns});
214 }
215 }
216 }
217
218
219 for (vector<hit_type>::const_iterator p = buffer.begin(); p != buffer.end(); ++p) {
220 for (vector<hit_type>::const_iterator q = buffer.begin(); p != q; ++q) {
221 h2r.Fill((double) combinatorics.getIndex(p->pmt,q->pmt),
222 JCombinatorics::getSign(p->pmt,q->pmt) * (q->t1 - p->t1));
223 }
224 }
225 }
226 }
227 STATUS(endl);
228
229
230 for (Int_t ix = 1; ix <= h2r.GetXaxis()->GetNbins(); ++ix) {
231 for (Int_t iy = 1; iy <= h2r.GetYaxis()->GetNbins(); ++iy) {
232 if (h2r.GetBinError(ix,iy) == 0.0) {
233 h2r.SetBinError(ix, iy, 1.0);
234 }
235 }
236 }
237
238 const double W = R / (WS * numberOfPrimaries);
239
240 h2r.Scale(W);
241
242 out << h2r;
243
244 out.Write();
245 out.Close();
246}
string outputFile
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define STATUS(A)
Definition JMessage.hh:63
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:72
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
#define MAKE_CSTRING(A)
Make C-string.
Definition JPrint.hh:72
Monte Carlo run header.
Definition JHead.hh:1236
JAANET::norma norma
Definition JHead.hh:1621
Detector data structure.
Definition JDetector.hh:96
Data structure for a composite optical module.
Definition JModule.hh:76
Template data structure for PMT I/O.
Data structure for PMT parameters.
double QE
relative quantum efficiency
double gainSpread
gain spread [unit]
double TTS_ns
transition time spread [ns]
JProperties getProperties(const JEquationParameters &equation=JPMTParameters::getEquationParameters())
Get properties of this class.
double PunderAmplified
probability of underamplified hit
bool slewing
time slewing of analogue signal
Utility class to parse parameter values.
General exception.
Definition JException.hh:24
int getID() const
Get identifier.
Definition JObjectID.hh:50
Utility class to parse command line options.
Definition JParser.hh:1698
General purpose class for object reading from a list of file names.
virtual bool hasNext() override
Check availability of next element.
counter_type getCounter() const
Get counter.
virtual const pointer_type & next() override
Get next element.
Auxiliary class to convert pair of indices to unique index and back.
static int getSign(const int first, const int second)
Sign of pair of indices.
Direct addressing of elements with unique identifiers.
Definition JRouter.hh:27
static const char *const _2R
Name extension for 2D rate measured.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
JCalibration getCalibration(const JCalibration &first, const JCalibration &second)
Get calibration to go from first to second calibration.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
bool putObject(TDirectory &dir, const TObject &object)
Write object to ROOT directory.
Head getHeader(const JMultipleFileScanner_t &file_list)
Get Monte Carlo header.
KM3NeT DAQ data structures and auxiliaries.
Definition DataQueue.cc:39
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition Evt.hh:21
std::vector< Hit > mc_hits
MC: list of MC truth hits.
Definition Evt.hh:48
The Head class reflects the header of Monte-Carlo event files, which consists of keys (also referred ...
Definition Head.hh:65
Detector file.
Definition JHead.hh:227
Transmission with position.
Definition JBillabong.cc:70
Auxiliary class to sort pairs of PMT addresses within optical module.
Data structure for PMT analogue signal.
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition JParser.hh:68
Auxiliary class for defining the range of iterations of objects.
Definition JLimit.hh:45
static counter_type max()
Get maximum counter value.
Definition JLimit.hh:128
Auxiliary class for ROOT I/O of application specific meta data.
Definition JMeta.hh:72