Jpp 19.3.0
the software that should make you happy
Loading...
Searching...
No Matches
JResolution.hh
Go to the documentation of this file.
1#ifndef __JASTRONOMY__JRESOLUTION__
2#define __JASTRONOMY__JRESOLUTION__
3
4#include <string>
5#include <istream>
6#include <ostream>
7#include <iomanip>
8#include <memory>
9#include <map>
10
11#include "TRandom3.h"
12#include "TROOT.h"
13#include "TFile.h"
14#include "TH1D.h"
15#include "TF1.h"
16
20
21#include "JLang/JTitle.hh"
22#include "JLang/JException.hh"
23
24
25/**
26 * \file
27 *
28 * Tools for simulation of detector resolution.
29 * \author mdejong
30 */
31namespace JASTRONOMY {}
32namespace JPP { using namespace JASTRONOMY; }
33
34namespace JASTRONOMY {
35
36 using JLANG::JTitle;
41
42
43 /**
44 * Inteface for detector resolution simulation.
45 * Note that all input angles are in degrees.
46 */
48 {
49 /**
50 * Virtual destructor.
51 */
52 virtual ~JResolution()
53 {}
54
55
56 /**
57 * Read detector resolution from input stream.
58 *
59 * \param in input stream
60 * \return input stream
61 */
62 virtual std::istream& read(std::istream& in) = 0;
63
64
65 /**
66 * Write detector resolution to output stream.
67 *
68 * \param out output stream
69 * \return output stream
70 */
71 virtual std::ostream& write(std::ostream& out) const = 0;
72
73
74 /**
75 * Get direction.
76 *
77 * \return direction
78 */
79 virtual JDirection3D get() const = 0;
80 };
81
82
83 /**
84 * Implementation of Gaussian detector resolution.
85 * Note that all input angles are in degrees.
86 */
88 public JResolution
89 {
90 /**
91 * Read detector resolution from input stream.
92 *
93 * \param in input stream
94 * \return input stream
95 */
96 virtual std::istream& read(std::istream& in) override
97 {
98 return in >> this->sigma_deg;
99 }
100
101
102 /**
103 * Write detector resolution to output stream.
104 *
105 * \param out output stream
106 * \return output stream
107 */
108 virtual std::ostream& write(std::ostream& out) const override
109 {
110 return out << this->sigma_deg;
111 }
112
113
114 /**
115 * Get direction.
116 *
117 * \return direction
118 */
119 JDirection3D get() const override
120 {
121 using namespace JPP;
122
123 JDirection3D u(0.0, 0.0, 1.0);
124
125 if (sigma_deg > 0.0) {
126
127 const double theta = gRandom->Gaus(0.0, getRadians(sigma_deg));
128 const double phi = gRandom->Uniform(-PI, +PI);
129
130 const JRotation3D R(JAngle3D(theta, phi));
131
132 u.rotate_back(R);
133 }
134
135 return u;
136 }
137
138
139 double sigma_deg;
140 };
141
142
143 /**
144 * Implementation of detector resolution based on a formula.
145 * Note that all input angles are in degrees.
146 */
148 public JResolution
149 {
150 /**
151 * Read detector resolution from input stream.
152 *
153 * \param in input stream
154 * \return input stream
155 */
156 virtual std::istream& read(std::istream& in) override
157 {
158 TString buffer;
159
160 if (buffer.ReadLine(in)) {
161
162 f1.reset(new TF1("f1", buffer));
163
164 if (!f1->IsValid()) {
165 THROW(JParseError, "Invalid function " << buffer);
166 }
167 }
168
169 return in;
170 }
171
172
173 /**
174 * Write detector resolution to output stream.
175 *
176 * \param out output stream
177 * \return output stream
178 */
179 virtual std::ostream& write(std::ostream& out) const override
180 {
181 if (f1)
182 return out << f1->GetExpFormula();
183 else
184 return out;
185 }
186
187
188 /**
189 * Get direction.
190 *
191 * \return direction
192 */
193 JDirection3D get() const override
194 {
195 using namespace JPP;
196
197 JDirection3D u(0.0, 0.0, 1.0);
198
199 if (f1) {
200
201 const double theta = getRadians(f1->GetRandom(gRandom));
202 const double phi = gRandom->Uniform(-PI, +PI);
203
204 const JRotation3D R(JAngle3D(theta, phi));
205
206 u.rotate_back(R);
207 }
208
209 return u;
210 }
211
212
213 std::shared_ptr<TF1> f1;
214 };
215
216
217 /**
218 * Implementation of detector resolution based on histogram.
219 * Note that all input angles are in degrees.
220 */
222 public JResolution
223 {
224 /**
225 * Virtual destructor.
226 */
228 {
229 if (in != NULL) {
230 in->Close();
231 }
232 }
233
234
235 /**
236 * Read detector resolution from input stream.
237 *
238 * \param in input stream
239 * \return input stream
240 */
241 virtual std::istream& read(std::istream& in) override
242 {
243 in >> this->filename
244 >> this->histname;
245
246 load();
247
248 return in;
249 }
250
251
252 /**
253 * Write detector resolution to output stream.
254 *
255 * \param out output stream
256 * \return output stream
257 */
258 virtual std::ostream& write(std::ostream& out) const override
259 {
260 return out << this->filename << ' '
261 << this->histname;
262 }
263
264
265 /**
266 * Get direction.
267 *
268 * \return direction
269 */
270 JDirection3D get() const override
271 {
272 using namespace JPP;
273
274 JDirection3D u(0.0, 0.0, 1.0);
275
276 if (h1 != NULL) {
277
278 const double theta = getRadians(h1->GetRandom(gRandom));
279 const double phi = gRandom->Uniform(-PI, +PI);
280
281 const JRotation3D R(JAngle3D(theta, phi));
282
283 u.rotate_back(R);
284 }
285
286 return u;
287 }
288
289
290 std::string filename;
291 std::string histname;
292
293 protected:
294 /**
295 * Load histogram from file.
296 */
297 void load()
298 {
299 in = TFile::Open(filename.c_str(), "exist");
300
301 if (in == NULL || !in->IsOpen()) {
302 THROW(JFileOpenException, "File: " << filename << " not opened.");
303 }
304
305 h1 = dynamic_cast<TH1D*>(in->Get(histname.c_str()));
306
307 if (h1 == NULL) {
308 THROW(JValueOutOfRange, "Histogram: " << histname << " not found.");
309 }
310 }
311
312 mutable TFile* in = NULL;
313 mutable TH1D* h1 = NULL;
314 };
315
316
317 /**
318 * Implementation of detector resolution based on histogram filled according <tt>log10(alpha)</tt>.
319 * Note that all input angles are in degrees.
320 */
322 public JResolutionTH1
323 {
324 /**
325 * Get direction.
326 *
327 * \return direction
328 */
329 JDirection3D get() const override
330 {
331 using namespace JPP;
332
333 JDirection3D u(0.0, 0.0, 1.0);
334
335 if (h1 != NULL) {
336
337 const double theta = getRadians(pow(10.0,h1->GetRandom(gRandom)));
338 const double phi = gRandom->Uniform(-PI, +PI);
339
340 const JRotation3D R(JAngle3D(theta, phi));
341
342 u.rotate_back(R);
343 }
344
345 return u;
346 }
347 };
348
349
350 /**
351 * Helper for detector resolution.
352 */
354 public JTitle,
355 public std::shared_ptr<JResolution>
356 {
357 /**
358 * Default constructor.
359 */
361 {
362 dictionary["Gauss"] .reset(new JResolutionGauss());
363 dictionary["Formula"] .reset(new JResolutionTF1());
364 dictionary["Histogram"] .reset(new JResolutionTH1());
365 dictionary["Logarithmic"].reset(new JResolutionLog());
366 }
367
368
369 /**
370 * Read detector resolution from input stream.
371 *
372 * \param in input stream
373 * \param object detector resolution
374 * \return input stream
375 */
376 friend inline std::istream& operator>>(std::istream& in, JResolution_t& object)
377 {
378 using namespace std;
379
380 string key;
381
382 if (in >> key) {
383
384 if (object.dictionary.count(key)) {
385
386 static_cast<std::shared_ptr<JResolution>&>(object) = object.dictionary[key];
387
388 object.setTitle(key);
389
390 object->read(in);
391
392 } else {
393
394 THROW(JParseError, "Invalid key " << key);
395 }
396 }
397
398 return in;
399 }
400
401
402 /**
403 * Write detector resolution to output stream.
404 *
405 * \param out output stream
406 * \param object detector resolution
407 * \return output stream
408 */
409 friend inline std::ostream& operator<<(std::ostream& out, const JResolution_t& object)
410 {
411 out << object.getTitle() << ' ';
412
413 if (static_cast<const std::shared_ptr<JResolution>&>(object)) {
414 object->write(out);
415 }
416
417 return out;
418 }
419
420 protected:
422 };
423}
424
425#endif
Exceptions.
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Data structure for angles in three dimensions.
Definition JAngle3D.hh:35
Data structure for direction in three dimensions.
JDirection3D & rotate_back(const JRotation3D &R)
Rotate back.
Exception for opening of file.
Exception for parsing value.
Auxiliary class for title.
Definition JTitle.hh:19
Exception for accessing a value in a collection that is outside of its range.
double getRadians(const double angle)
Convert angle to radians.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Implementation of Gaussian detector resolution.
virtual std::istream & read(std::istream &in) override
Read detector resolution from input stream.
JDirection3D get() const override
Get direction.
virtual std::ostream & write(std::ostream &out) const override
Write detector resolution to output stream.
Implementation of detector resolution based on histogram filled according <tt>log10(alpha).
JDirection3D get() const override
Get direction.
Implementation of detector resolution based on a formula.
virtual std::istream & read(std::istream &in) override
Read detector resolution from input stream.
JDirection3D get() const override
Get direction.
std::shared_ptr< TF1 > f1
virtual std::ostream & write(std::ostream &out) const override
Write detector resolution to output stream.
Implementation of detector resolution based on histogram.
virtual std::istream & read(std::istream &in) override
Read detector resolution from input stream.
virtual ~JResolutionTH1()
Virtual destructor.
void load()
Load histogram from file.
JDirection3D get() const override
Get direction.
virtual std::ostream & write(std::ostream &out) const override
Write detector resolution to output stream.
Helper for detector resolution.
friend std::ostream & operator<<(std::ostream &out, const JResolution_t &object)
Write detector resolution to output stream.
JResolution_t()
Default constructor.
friend std::istream & operator>>(std::istream &in, JResolution_t &object)
Read detector resolution from input stream.
std::map< std::string, std::shared_ptr< JResolution > > dictionary
Inteface for detector resolution simulation.
virtual JDirection3D get() const =0
Get direction.
virtual ~JResolution()
Virtual destructor.
virtual std::istream & read(std::istream &in)=0
Read detector resolution from input stream.
virtual std::ostream & write(std::ostream &out) const =0
Write detector resolution to output stream.