Jpp test-rotations-old-533-g2bdbdb559
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
25
26
27/**
28 * \file
29 *
30 * Tools for simulation of detector resolution.
31 * \author mdejong
32 */
33namespace JASTRONOMY {}
34namespace JPP { using namespace JASTRONOMY; }
35
36namespace JASTRONOMY {
37
38 using JLANG::JTitle;
43
44
45 /**
46 * Interface for detector resolution simulation.
47 * Note that all input angles are in degrees.
48 */
50 {
51 /**
52 * Virtual destructor.
53 */
54 virtual ~JResolution()
55 {}
56
57
58 /**
59 * Read detector resolution from input stream.
60 *
61 * \param in input stream
62 * \return input stream
63 */
64 virtual std::istream& read(std::istream& in) = 0;
65
66
67 /**
68 * Write detector resolution to output stream.
69 *
70 * \param out output stream
71 * \return output stream
72 */
73 virtual std::ostream& write(std::ostream& out) const = 0;
74
75
76 /**
77 * Get direction.
78 *
79 * \return direction
80 */
81 virtual JDirection3D get() const = 0;
82 };
83
84
85 /**
86 * Implementation of Gaussian detector resolution.
87 * Note that all input angles are in degrees.
88 */
90 public JResolution
91 {
92 /**
93 * Read detector resolution from input stream.
94 *
95 * \param in input stream
96 * \return input stream
97 */
98 virtual std::istream& read(std::istream& in) override
99 {
100 return in >> this->sigma_deg;
101 }
102
103
104 /**
105 * Write detector resolution to output stream.
106 *
107 * \param out output stream
108 * \return output stream
109 */
110 virtual std::ostream& write(std::ostream& out) const override
111 {
112 return out << this->sigma_deg;
113 }
114
115
116 /**
117 * Get direction.
118 *
119 * \return direction
120 */
121 JDirection3D get() const override
122 {
123 using namespace JPP;
124
125 JDirection3D u(0.0, 0.0, 1.0);
126
127 if (sigma_deg > 0.0) {
128
129 const double theta = gRandom->Gaus(0.0, getRadians(sigma_deg));
130 const double phi = gRandom->Uniform(-PI, +PI);
131
132 const JRotation3D R(JAngle3D(theta, phi));
133
134 u.rotate_back(R);
135 }
136
137 return u;
138 }
139
140
141 double sigma_deg;
142 };
143
144
145 /**
146 * Implementation of detector resolution based on a formula.
147 * Note that all input angles are in degrees.
148 */
150 public JResolution
151 {
152 /**
153 * Read detector resolution from input stream.
154 *
155 * \param in input stream
156 * \return input stream
157 */
158 virtual std::istream& read(std::istream& in) override
159 {
160 TString buffer;
161
162 if (buffer.ReadLine(in)) {
163
164 f1.reset(new TF1("f1", buffer));
165
166 if (!f1->IsValid()) {
167 THROW(JParseError, "Invalid function " << buffer);
168 }
169 }
170
171 return in;
172 }
173
174
175 /**
176 * Write detector resolution to output stream.
177 *
178 * \param out output stream
179 * \return output stream
180 */
181 virtual std::ostream& write(std::ostream& out) const override
182 {
183 if (f1)
184 return out << f1->GetExpFormula();
185 else
186 return out;
187 }
188
189
190 /**
191 * Get direction.
192 *
193 * \return direction
194 */
195 JDirection3D get() const override
196 {
197 using namespace JPP;
198
199 JDirection3D u(0.0, 0.0, 1.0);
200
201 if (f1) {
202
203 const double theta = getRadians(f1->GetRandom(gRandom));
204 const double phi = gRandom->Uniform(-PI, +PI);
205
206 const JRotation3D R(JAngle3D(theta, phi));
207
208 u.rotate_back(R);
209 }
210
211 return u;
212 }
213
214
215 std::shared_ptr<TF1> f1;
216 };
217
218
219 /**
220 * Implementation of detector resolution based on histogram.
221 * Note that all input angles are in degrees.
222 */
224 public JResolution
225 {
226 /**
227 * Virtual destructor.
228 */
230 {
231 if (in != NULL) {
232 in->Close();
233 }
234 }
235
236
237 /**
238 * Read detector resolution from input stream.
239 *
240 * \param in input stream
241 * \return input stream
242 */
243 virtual std::istream& read(std::istream& in) override
244 {
245 in >> this->filename
246 >> this->histname;
247
248 load();
249
250 return in;
251 }
252
253
254 /**
255 * Write detector resolution to output stream.
256 *
257 * \param out output stream
258 * \return output stream
259 */
260 virtual std::ostream& write(std::ostream& out) const override
261 {
262 return out << this->filename << ' '
263 << this->histname;
264 }
265
266
267 /**
268 * Get direction.
269 *
270 * \return direction
271 */
272 JDirection3D get() const override
273 {
274 using namespace JPP;
275
276 JDirection3D u(0.0, 0.0, 1.0);
277
278 if (h1 != NULL) {
279
280 const double theta = getRadians(h1->GetRandom(gRandom));
281 const double phi = gRandom->Uniform(-PI, +PI);
282
283 const JRotation3D R(JAngle3D(theta, phi));
284
285 u.rotate_back(R);
286 }
287
288 return u;
289 }
290
291
292 std::string filename;
293 std::string histname;
294
295 protected:
296 /**
297 * Load histogram from file.
298 */
299 void load()
300 {
301 in = TFile::Open(filename.c_str(), "exist");
302
303 if (in == NULL || !in->IsOpen()) {
304 THROW(JFileOpenException, "File: " << filename << " not opened.");
305 }
306
307 h1 = dynamic_cast<TH1D*>(in->Get(histname.c_str()));
308
309 if (h1 == NULL) {
310 THROW(JValueOutOfRange, "Histogram: " << histname << " not found.");
311 }
312 }
313
314 mutable TFile* in = NULL;
315 mutable TH1D* h1 = NULL;
316 };
317
318
319 /**
320 * Implementation of detector resolution based on histogram filled according <tt>log10(alpha)</tt>.
321 * Note that all input angles are in degrees.
322 */
324 public JResolutionTH1
325 {
326 /**
327 * Get direction.
328 *
329 * \return direction
330 */
331 JDirection3D get() const override
332 {
333 using namespace JPP;
334
335 JDirection3D u(0.0, 0.0, 1.0);
336
337 if (h1 != NULL) {
338
339 const double theta = getRadians(pow(10.0,h1->GetRandom(gRandom)));
340 const double phi = gRandom->Uniform(-PI, +PI);
341
342 const JRotation3D R(JAngle3D(theta, phi));
343
344 u.rotate_back(R);
345 }
346
347 return u;
348 }
349 };
350
351
352 /**
353 * Type definition of generic resolution.
354 */
355 typedef std::shared_ptr<JResolution> resolution_type;
356
357
358 /**
359 * Helper for detector resolution I/O.
360 */
362 public JTitle,
363 public resolution_type
364 {
366
367 /**
368 * Default constructor.
369 */
371 {
372 dictionary["Gauss"] .reset(new JResolutionGauss());
373 dictionary["Formula"] .reset(new JResolutionTF1());
374 dictionary["Histogram"] .reset(new JResolutionTH1());
375 dictionary["Logarithmic"].reset(new JResolutionLog());
376 }
377
378
379 /**
380 * Get dictionary.
381 *
382 * \return dictionary
383 */
385 {
386 return dictionary;
387 }
388
389
390 /**
391 * Read detector resolution from input stream.
392 *
393 * \param in input stream
394 * \param object detector resolution
395 * \return input stream
396 */
397 friend inline std::istream& operator>>(std::istream& in, JResolution_t& object)
398 {
399 using namespace std;
400
401 string key;
402
403 if (in >> key) {
404
405 if (object.dictionary.count(key)) {
406
407 static_cast<shared_ptr<JResolution>&>(object) = object.dictionary[key];
408
409 object.setTitle(key);
410
411 object->read(in);
412
413 } else {
414
415 THROW(JParseError, "Invalid key " << key);
416 }
417 }
418
419 return in;
420 }
421
422
423 /**
424 * Write detector resolution to output stream.
425 *
426 * \param out output stream
427 * \param object detector resolution
428 * \return output stream
429 */
430 friend inline std::ostream& operator<<(std::ostream& out, const JResolution_t& object)
431 {
432 out << object.getTitle() << ' ';
433
434 if (static_cast<const std::shared_ptr<JResolution>&>(object)) {
435 object->write(out);
436 }
437
438 return out;
439 }
440
441 protected:
443 };
444}
445
446#endif
Interface methods for SLALIB and auxiliary classes and methods for astronomy.
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.
std::shared_ptr< JResolution > resolution_type
Type definition of generic resolution.
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 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 I/O.
dictionary_type dictionary
std::map< std::string, resolution_type > dictionary_type
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.
const dictionary_type & getDictionary() const
Get dictionary.
Interface 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.