Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
DoubleAc.h
Go to the documentation of this file.
1#ifndef DOUBLEAC_H
2#define DOUBLEAC_H
3
4/*
5Copyright (c) 2001 Igor B. Smirnov
6
7The file can be used, copied, modified, and distributed
8according to the terms of GNU Lesser General Public License version 2.1
9as published by the Free Software Foundation,
10and provided that the above copyright notice, this permission notice,
11and notices about any modifications of the original text
12appear in all copies and in supporting documentation.
13The file is provided "as is" without express or implied warranty.
14*/
15
16#ifdef VISUAL_STUDIO
17#define _USE_MATH_DEFINES
18/* Define _USE_MATH_DEFINES before including math.h to expose these macro
19 * definitions for common math constants. These are placed under an #ifdef
20 * since these commonly-defined names are not part of the C/C++ standards.
21 */
22#endif
23#include <cmath>
24#include <climits>
25#include <cfloat>
27
28namespace Heed {
29
30#define DEF_DBL_PREC 1.0e-15 // rounding precision for double
31const double one_plus_def_dbl_prec = double(1.0) + DEF_DBL_PREC;
32const double one_minus_def_dbl_prec = double(1.0) - DEF_DBL_PREC;
33
34#define DEF_FLT_PREC 1.0e-7 // rounding precision for float
35const double one_plus_def_flt_prec = double(1.0) + DEF_FLT_PREC;
36const double one_minus_def_flt_prec = double(1.0) - DEF_FLT_PREC;
37
38/// "Double with accuracy".
39/// Really this is original implementation of interval computations.
40/// This is the way of watching for the numerical errors
41/// by establishing lower and upper limit of each numerical value.
42/// The central the most probable value is also watched for.
43
44class DoubleAc {
45 double d; //< the main value
46 double di; //< the left limit di <= d
47 double da; //< the right limit da >= d
48 public:
49 inline DoubleAc(void) {
50 di = 0.0;
51 d = 0.0;
52 da = 0.0;
53 }
54 inline DoubleAc(const DoubleAc& f) {
55 di = f.di;
56 d = f.d;
57 da = f.da;
58 }
59 inline DoubleAc(double f);
60 inline DoubleAc(float f);
61 inline DoubleAc(long f);
62 inline DoubleAc(int f);
63 DoubleAc(double f, double ffmin, double ffmax);
64 DoubleAc(double f, double relative_prec);
65 inline DoubleAc& operator=(const DoubleAc& f) {
66 d = f.d;
67 di = f.di;
68 da = f.da;
69 return *this;
70 }
71 inline DoubleAc& operator=(double f);
72 inline DoubleAc& operator=(float f);
73 inline DoubleAc& operator=(long f);
74 inline DoubleAc& operator=(int f);
75 inline operator double(void) const { return d; }
76 inline double get(void) const { return d; }
77 inline double get_low_limit(void) const { return di; }
78 inline double get_left_limit(void) const { return di; }
79 inline double left_limit(void) const { return di; }
80 inline double get_min_limit(void) const { return di; }
81 inline double get_high_limit(void) const { return da; }
82 inline double get_right_limit(void) const { return da; }
83 inline double right_limit(void) const { return da; }
84 inline double get_max_limit(void) const { return da; }
85 inline double get_accuracy(void) const { return 0.5 * (da - di); }
86 inline DoubleAc& operator+=(const DoubleAc& f);
87 inline DoubleAc& operator+=(double f);
88 inline DoubleAc& operator+=(float f);
89 inline DoubleAc& operator+=(long f);
90 inline DoubleAc& operator+=(int f);
91
92 inline DoubleAc& operator-=(const DoubleAc& f);
93 inline DoubleAc& operator-=(double f);
94 inline DoubleAc& operator-=(float f);
95 inline DoubleAc& operator-=(long f);
96 inline DoubleAc& operator-=(int f);
97
98 friend inline void change_sign(DoubleAc& f);
99
101 inline DoubleAc& operator*=(double f);
102 inline DoubleAc& operator*=(float f);
103 inline DoubleAc& operator*=(long f);
104 inline DoubleAc& operator*=(int f);
105
107 inline DoubleAc& operator/=(double f);
108 inline DoubleAc& operator/=(float f);
109 inline DoubleAc& operator/=(long f);
110 inline DoubleAc& operator/=(int f);
111
112 void print(std::ostream& file, int l = 1) const;
113 // if(l <= 0) return without print
114 // if i == 1, 2, or 3, print without passing to the next line at the end
115 // if i == 4, 5, or 6, the same print but with passing to the next line
116 // (sending file<<'\n';)
117 // 1 - output to file only central value by ordinary operator file<<d;
118 // 2 - output central value and both limits in brackets.
119 // 3 - output central value, then set e.precision(2),
120 // output differences in brackets d-di and da-d (both positive),
121 // and restore initial precision.
122 // 4 - output to file only central value by ordinary operator file<<d;
123 // and passing to the next line.
124 // 5 - the same as 3 but passing to the next line at the end.
125 // 6 - output all values with precision 16 and with width 20,
126 // restore initial precision, and passing to the next line.
127 // The most convenient options are 3 or 5.
128};
129
130inline DoubleAc::DoubleAc(double f) {
131 // assumes pure numerical uncertanty which is defined via macro above
132 // unless f == 0. In this case the precision is assumed absolute.
133 d = f;
134 if (f == 0.0) {
135 di = 0.0;
136 da = 0.0;
137 } else {
138 if (f > 0.0) {
139 if (f < DBL_MAX / one_plus_def_dbl_prec)
140 da = f * one_plus_def_dbl_prec;
141 else
142 da = f;
143 if (f > DBL_MIN * one_plus_def_dbl_prec)
144 di = f / one_plus_def_dbl_prec;
145 else
146 di = f;
147 } else {
148 if (-f < DBL_MAX / one_plus_def_dbl_prec)
149 di = f * one_plus_def_dbl_prec;
150 else
151 di = f;
152 if (-f > DBL_MIN * one_plus_def_dbl_prec)
153 da = f / one_plus_def_dbl_prec;
154 else
155 da = f;
156 }
157 }
158}
159
160inline DoubleAc& DoubleAc::operator=(double f) {
161 d = f;
162 if (f == 0.0) {
163 di = 0.0;
164 da = 0.0;
165 } else {
166 if (f > 0.0) {
167 if (f < DBL_MAX / one_plus_def_dbl_prec)
168 da = f * one_plus_def_dbl_prec;
169 else
170 da = f;
171 if (f > DBL_MIN * one_plus_def_dbl_prec)
172 di = f / one_plus_def_dbl_prec;
173 else
174 di = f;
175 } else {
176 if (-f < DBL_MAX / one_plus_def_dbl_prec)
177 di = f * one_plus_def_dbl_prec;
178 else
179 di = f;
180 if (-f > DBL_MIN * one_plus_def_dbl_prec)
181 da = f / one_plus_def_dbl_prec;
182 else
183 da = f;
184 }
185 }
186 return *this;
187}
188
189inline DoubleAc::DoubleAc(float f) {
190 d = f;
191 if (f == 0.0) {
192 di = 0.0;
193 da = 0.0;
194 } else if (f > 0.0) {
195 di = f * one_minus_def_flt_prec;
196 da = f * one_plus_def_flt_prec;
197 } else {
198 da = f * one_minus_def_flt_prec;
199 di = f * one_plus_def_flt_prec;
200 }
201}
202
204 d = f;
205 if (f == 0.0) {
206 di = 0.0;
207 da = 0.0;
208 } else if (f > 0.0) {
209 di = f * one_minus_def_flt_prec;
210 da = f * one_plus_def_flt_prec;
211 } else {
212 da = f * one_minus_def_flt_prec;
213 di = f * one_plus_def_flt_prec;
214 }
215 return *this;
216}
217
218inline DoubleAc::DoubleAc(long f) {
219 d = f;
220 di = f;
221 da = f;
222}
223
225 d = f;
226 di = f;
227 da = f;
228 return *this;
229}
230
231inline DoubleAc::DoubleAc(int f) {
232 d = f;
233 di = f;
234 da = f;
235}
237 d = f;
238 di = f;
239 da = f;
240 return *this;
241}
242
244 di += f.di;
245 d += f.d;
246 da += f.da;
247 return *this;
248}
250 *this += DoubleAc(f);
251 return *this;
252}
253//{ di+=f; d+=f; da+=f; return *this; }
255 *this += DoubleAc(f);
256 return *this;
257}
258//{ di+=f; d+=f; da+=f; return *this; }
260 di += f;
261 d += f;
262 da += f;
263 return *this;
264}
266 di += f;
267 d += f;
268 da += f;
269 return *this;
270}
271
273 di -= f.da;
274 d -= f.d;
275 da -= f.di;
276 return *this;
277}
279 *this -= DoubleAc(f);
280 return *this;
281}
282//{ di-=f; d-=f; da-=f; return *this; }
284 *this -= DoubleAc(f);
285 return *this;
286}
287//{ di-=f; d-=f; da-=f; return *this; }
289 di -= f;
290 d -= f;
291 da -= f;
292 return *this;
293}
295 di -= f;
296 d -= f;
297 da -= f;
298 return *this;
299}
300
302 d *= f;
303 if (f >= 0.0) {
304 di *= f;
305 da *= f;
306 } else {
307 double ti = da * f;
308 da = di * f;
309 di = ti;
310 }
311 return *this;
312}
314 d *= f;
315 if (f >= 0.0) {
316 di *= f;
317 da *= f;
318 } else {
319 double ti = da * f;
320 da = di * f;
321 di = ti;
322 }
323 return *this;
324}
326 d *= f;
327 if (f >= 0.0) {
328 di *= f;
329 da *= f;
330 } else {
331 double ti = da * f;
332 da = di * f;
333 di = ti;
334 }
335 return *this;
336}
338 d *= f;
339 if (f >= 0.0) {
340 di *= f;
341 da *= f;
342 } else {
343 double ti = da * f;
344 da = di * f;
345 di = ti;
346 }
347 return *this;
348}
349
351 if (f == 0.0) {
352 mcerr << "inline DoubleAc& DoubleAc::operator/=(double f):\n"
353 << "f = 0.0\n";
354 spexit(mcerr);
355 }
356 d /= f;
357 if (f >= 0.0) {
358 di /= f;
359 da /= f;
360 } else {
361 double ti = da / f;
362 da = di / f;
363 di = ti;
364 }
365 return *this;
366}
368 if (f == 0.0) {
369 mcerr << "inline DoubleAc& DoubleAc::operator/=(float f):\n"
370 << "f = 0.0\n";
371 spexit(mcerr);
372 }
373 d /= f;
374 if (f >= 0.0) {
375 di /= f;
376 da /= f;
377 } else {
378 double ti = da / f;
379 da = di / f;
380 di = ti;
381 }
382 return *this;
383}
385 if (f == 0) {
386 mcerr << "inline DoubleAc& DoubleAc::operator/=(long f):\n"
387 << "f = 0\n";
388 spexit(mcerr);
389 }
390 d /= f;
391 if (f >= 0.0) {
392 di /= f;
393 da /= f;
394 } else {
395 double ti = da / f;
396 da = di / f;
397 di = ti;
398 }
399 return *this;
400}
402 if (f == 0) {
403 mcerr << "inline DoubleAc& DoubleAc::operator/=(int f):\n"
404 << "f = 0\n";
405 spexit(mcerr);
406 }
407 d /= f;
408 if (f >= 0.0) {
409 di /= f;
410 da /= f;
411 } else {
412 double ti = da / f;
413 da = di / f;
414 di = ti;
415 }
416 return *this;
417}
418
419inline DoubleAc operator-(const DoubleAc& f) {
420 DoubleAc t(-f.get(), -f.get_right_limit(), -f.get_left_limit());
421 return t;
422}
423
424inline void change_sign(DoubleAc& f) {
425 f.d = -f.d;
426 double temp = f.di;
427 f.di = -f.da;
428 f.da = -temp;
429}
430
431inline DoubleAc operator+(const DoubleAc& f1, const DoubleAc& f2) {
432 DoubleAc t = f1;
433 t += f2;
434 return t;
435}
436inline DoubleAc operator+(const DoubleAc& f1, double f2) {
437 DoubleAc t = f1;
438 t += f2;
439 return t;
440}
441inline DoubleAc operator+(double f1, const DoubleAc& f2) {
442 DoubleAc t = f2;
443 t += f1;
444 return t;
445}
446inline DoubleAc operator+(const DoubleAc& f1, float f2) {
447 DoubleAc t = f1;
448 t += f2;
449 return t;
450}
451inline DoubleAc operator+(float f1, const DoubleAc& f2) {
452 DoubleAc t = f2;
453 t += f1;
454 return t;
455}
456inline DoubleAc operator+(const DoubleAc& f1, long f2) {
457 DoubleAc t = f1;
458 t += f2;
459 return t;
460}
461inline DoubleAc operator+(long f1, const DoubleAc& f2) {
462 DoubleAc t = f2;
463 t += f1;
464 return t;
465}
466inline DoubleAc operator+(const DoubleAc& f1, int f2) {
467 DoubleAc t = f1;
468 t += f2;
469 return t;
470}
471inline DoubleAc operator+(int f1, const DoubleAc& f2) {
472 DoubleAc t = f2;
473 t += f1;
474 return t;
475}
476
477inline DoubleAc operator-(const DoubleAc& f1, const DoubleAc& f2) {
478 DoubleAc t = f1;
479 t -= f2;
480 return t;
481}
482inline DoubleAc operator-(const DoubleAc& f1, double f2) {
483 DoubleAc t = f1;
484 t -= f2;
485 return t;
486}
487inline DoubleAc operator-(double f1, const DoubleAc& f2) {
488 DoubleAc t = -f2;
489 t += f1;
490 return t;
491}
492inline DoubleAc operator-(const DoubleAc& f1, float f2) {
493 DoubleAc t = f1;
494 t -= f2;
495 return t;
496}
497inline DoubleAc operator-(float f1, const DoubleAc& f2) {
498 DoubleAc t = -f2;
499 t += f1;
500 return t;
501}
502inline DoubleAc operator-(const DoubleAc& f1, long f2) {
503 DoubleAc t = f1;
504 t -= f2;
505 return t;
506}
507inline DoubleAc operator-(long f1, const DoubleAc& f2) {
508 DoubleAc t = -f2;
509 t += f1;
510 return t;
511}
512inline DoubleAc operator-(const DoubleAc& f1, int f2) {
513 DoubleAc t = f1;
514 t -= f2;
515 return t;
516}
517inline DoubleAc operator-(int f1, const DoubleAc& f2) {
518 DoubleAc t = -f2;
519 t += f1;
520 return t;
521}
522
523inline DoubleAc operator*(const DoubleAc& f1, const DoubleAc& f2) {
524 DoubleAc t = f1;
525 t *= f2;
526 return t;
527}
528inline DoubleAc operator*(const DoubleAc& f1, double f2) {
529 DoubleAc t = f1;
530 t *= f2;
531 return t;
532}
533inline DoubleAc operator*(double f1, const DoubleAc& f2) {
534 DoubleAc t = f2;
535 t *= f1;
536 return t;
537}
538inline DoubleAc operator*(const DoubleAc& f1, float f2) {
539 DoubleAc t = f1;
540 t *= f2;
541 return t;
542}
543inline DoubleAc operator*(float f1, const DoubleAc& f2) {
544 DoubleAc t = f2;
545 t *= f1;
546 return t;
547}
548inline DoubleAc operator*(const DoubleAc& f1, long f2) {
549 DoubleAc t = f1;
550 t *= f2;
551 return t;
552}
553inline DoubleAc operator*(long f1, const DoubleAc& f2) {
554 DoubleAc t = f2;
555 t *= f1;
556 return t;
557}
558inline DoubleAc operator*(const DoubleAc& f1, int f2) {
559 DoubleAc t = f1;
560 t *= f2;
561 return t;
562}
563inline DoubleAc operator*(int f1, const DoubleAc& f2) {
564 DoubleAc t = f2;
565 t *= f1;
566 return t;
567}
568
569inline DoubleAc operator/(const DoubleAc& f1, const DoubleAc& f2) {
570 DoubleAc t = f1;
571 t /= f2;
572 return t;
573}
574inline DoubleAc operator/(const DoubleAc& f1, double f2) {
575 DoubleAc t = f1;
576 t /= f2;
577 return t;
578}
579inline DoubleAc operator/(double f1, const DoubleAc& f2) {
580 DoubleAc t = f1;
581 t /= f2;
582 return t;
583}
584inline DoubleAc operator/(const DoubleAc& f1, float f2) {
585 DoubleAc t = f1;
586 t /= f2;
587 return t;
588}
589inline DoubleAc operator/(float f1, const DoubleAc& f2) {
590 DoubleAc t = f1;
591 t /= f2;
592 return t;
593}
594inline DoubleAc operator/(const DoubleAc& f1, long f2) {
595 DoubleAc t = f1;
596 t /= f2;
597 return t;
598}
599inline DoubleAc operator/(long f1, const DoubleAc& f2) {
600 DoubleAc t = f1;
601 t /= f2;
602 return t;
603}
604inline DoubleAc operator/(const DoubleAc& f1, int f2) {
605 DoubleAc t = f1;
606 t /= f2;
607 return t;
608}
609inline DoubleAc operator/(int f1, const DoubleAc& f2) {
610 DoubleAc t = f1;
611 t /= f2;
612 return t;
613}
614
615inline DoubleAc fabs(const DoubleAc& f) {
616 if (f.left_limit() >= 0)
617 return f;
618 else if (f.right_limit() > 0) {
619 return DoubleAc(fabs(f.get()), 0.0,
620 std::max(f.right_limit(), -f.left_limit()));
621 } else {
622 return DoubleAc(-f.get(), -f.right_limit(), -f.left_limit());
623 }
624}
625
626inline DoubleAc find_min(const DoubleAc& a, const DoubleAc& b) {
627 return (a.get() < b.get() ? a : b);
628}
629inline DoubleAc find_min(const DoubleAc& a, double b) {
630 return (a.get() < b ? a : DoubleAc(b));
631}
632inline DoubleAc find_min(double a, const DoubleAc& b) {
633 return (a < b.get() ? DoubleAc(a) : b);
634}
635inline DoubleAc find_min(const DoubleAc& a, float b) {
636 return (a.get() < b ? a : DoubleAc(b));
637}
638inline DoubleAc find_min(float a, const DoubleAc& b) {
639 return (a < b.get() ? DoubleAc(a) : b);
640}
641inline DoubleAc find_min(const DoubleAc& a, long b) {
642 return (a.get() < b ? a : DoubleAc(b));
643}
644inline DoubleAc find_min(long a, const DoubleAc& b) {
645 return (a < b.get() ? DoubleAc(a) : b);
646}
647inline DoubleAc find_min(const DoubleAc& a, int b) {
648 return (a.get() < b ? a : DoubleAc(b));
649}
650inline DoubleAc find_min(int a, const DoubleAc& b) {
651 return (a < b.get() ? DoubleAc(a) : b);
652}
653
654inline DoubleAc find_max(const DoubleAc& a, const DoubleAc& b) {
655 return (a.get() > b.get() ? a : b);
656}
657inline DoubleAc find_max(const DoubleAc& a, double b) {
658 return (a.get() > b ? a : DoubleAc(b));
659}
660inline DoubleAc find_max(double a, const DoubleAc& b) {
661 return (a > b.get() ? DoubleAc(a) : b);
662}
663inline DoubleAc find_max(const DoubleAc& a, float b) {
664 return (a.get() > b ? a : DoubleAc(b));
665}
666inline DoubleAc find_max(float a, const DoubleAc& b) {
667 return (a > b.get() ? DoubleAc(a) : b);
668}
669inline DoubleAc find_max(const DoubleAc& a, long b) {
670 return (a.get() > b ? a : DoubleAc(b));
671}
672inline DoubleAc find_max(long a, const DoubleAc& b) {
673 return (a > b.get() ? DoubleAc(a) : b);
674}
675inline DoubleAc find_max(const DoubleAc& a, int b) {
676 return (a.get() > b ? a : DoubleAc(b));
677}
678inline DoubleAc find_max(int a, const DoubleAc& b) {
679 return (a > b.get() ? DoubleAc(a) : b);
680}
681
682DoubleAc sqrt(const DoubleAc& f);
683
684DoubleAc square(const DoubleAc& f);
685
686DoubleAc pow(const DoubleAc& f, double p);
687
688DoubleAc exp(const DoubleAc& f);
689
690DoubleAc sin(const DoubleAc& f);
691
692DoubleAc cos(const DoubleAc& f);
693
694DoubleAc asin(const DoubleAc& f);
695
696DoubleAc acos(const DoubleAc& f);
697
698std::ostream& operator<<(std::ostream& file, const DoubleAc& f);
699// Calls f.print(file, 1) (prints only the central value)
700
701#define Iprintda(file, name) \
702 file << indn << #name << "=" << noindent; \
703 name.print(file, 3); \
704 file << yesindent;
705
706#define Iprintdan(file, name) \
707 file << indn << #name << "=" << noindent; \
708 name.print(file, 3); \
709 file << yesindent << '\n';
710
711}
712
713#endif
#define DEF_FLT_PREC
Definition: DoubleAc.h:34
#define DEF_DBL_PREC
Definition: DoubleAc.h:30
#define spexit(stream)
Definition: FunNameStack.h:256
DoubleAc & operator=(const DoubleAc &f)
Definition: DoubleAc.h:65
double get_right_limit(void) const
Definition: DoubleAc.h:82
DoubleAc & operator/=(DoubleAc f)
Definition: DoubleAc.cpp:190
double get_low_limit(void) const
Definition: DoubleAc.h:77
DoubleAc(void)
Definition: DoubleAc.h:49
void print(std::ostream &file, int l=1) const
Definition: DoubleAc.cpp:510
double get(void) const
Definition: DoubleAc.h:76
friend void change_sign(DoubleAc &f)
Definition: DoubleAc.h:424
double get_left_limit(void) const
Definition: DoubleAc.h:78
DoubleAc & operator-=(const DoubleAc &f)
Definition: DoubleAc.h:272
double left_limit(void) const
Definition: DoubleAc.h:79
DoubleAc(const DoubleAc &f)
Definition: DoubleAc.h:54
double right_limit(void) const
Definition: DoubleAc.h:83
DoubleAc & operator+=(const DoubleAc &f)
Definition: DoubleAc.h:243
double get_min_limit(void) const
Definition: DoubleAc.h:80
DoubleAc & operator*=(DoubleAc f)
Definition: DoubleAc.cpp:41
double get_high_limit(void) const
Definition: DoubleAc.h:81
double get_accuracy(void) const
Definition: DoubleAc.h:85
double get_max_limit(void) const
Definition: DoubleAc.h:84
Definition: BGMesh.cpp:6
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
DoubleAc operator+(const DoubleAc &f1, const DoubleAc &f2)
Definition: DoubleAc.h:431
std::ostream & operator<<(std::ostream &file, const BGMesh &bgm)
Definition: BGMesh.cpp:37
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:337
DoubleAc operator-(const DoubleAc &f)
Definition: DoubleAc.h:419
DoubleAc operator/(const DoubleAc &f1, const DoubleAc &f2)
Definition: DoubleAc.h:569
DoubleAc exp(const DoubleAc &f)
Definition: DoubleAc.cpp:377
DoubleAc find_min(const DoubleAc &a, const DoubleAc &b)
Definition: DoubleAc.h:626
DoubleAc square(const DoubleAc &f)
Definition: DoubleAc.cpp:325
void change_sign(DoubleAc &f)
Definition: DoubleAc.h:424
DoubleAc acos(const DoubleAc &f)
Definition: DoubleAc.cpp:490
DoubleAc operator*(const DoubleAc &f1, const DoubleAc &f2)
Definition: DoubleAc.h:523
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384
DoubleAc find_max(const DoubleAc &a, const DoubleAc &b)
Definition: DoubleAc.h:654
DoubleAc asin(const DoubleAc &f)
Definition: DoubleAc.cpp:470
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314
const double one_plus_def_flt_prec
Definition: DoubleAc.h:35
const double one_minus_def_dbl_prec
Definition: DoubleAc.h:32
const double one_plus_def_dbl_prec
Definition: DoubleAc.h:31
const double one_minus_def_flt_prec
Definition: DoubleAc.h:36
#define mcerr
Definition: prstream.h:128