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