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.cpp
Go to the documentation of this file.
1/*
2Copyright (c) 2001 Igor B. Smirnov
3
4The file can be used, copied, modified, and distributed
5according to the terms of GNU Lesser General Public License version 2.1
6as published by the Free Software Foundation,
7and provided that the above copyright notice, this permission notice,
8and notices about any modifications of the original text
9appear in all copies and in supporting documentation.
10The file is provided "as is" without express or implied warranty.
11*/
12#include <iomanip>
13#include "wcpplib/math/minmax.h"
15
16namespace Heed {
17
18DoubleAc::DoubleAc(double f, double ffmin, double ffmax) {
19 mfunname("DoubleAc::DoubleAc(double f, double ffmin, double ffmax)");
20 check_econd12a(f, <, ffmin, "f - ffmin=" << f - ffmin << '\n', mcerr);
21 check_econd12a(f, >, ffmax, "f - ffmax=" << f - ffmax << '\n', mcerr);
22 di = ffmin;
23 d = f;
24 da = ffmax;
25}
26
27DoubleAc::DoubleAc(double f, double relative_prec) {
28 mfunname("DoubleAc::DoubleAc(double f, double relative_prec)");
29 check_econd11(relative_prec, < 0, mcerr);
30 check_econd11(relative_prec, >= 1, mcerr);
31 d = f;
32 if (f >= 0) {
33 da = f * (1.0 + relative_prec);
34 di = f / (1.0 + relative_prec);
35 } else {
36 di = f * (1.0 + relative_prec);
37 da = f / (1.0 + relative_prec);
38 }
39}
40
42 mfunnamep("DoubleAc& DoubleAc::operator*=(const DoubleAc& f)");
43#ifdef DEBUG_OPERATOR_MULT_PRINT
44 mcout.precision(17);
45 mcout << "d=" << d << '\n';
46 mcout << "di=" << di << '\n';
47 mcout << "da=" << da << '\n';
48 mcout << "f.d=" << f.d << '\n';
49 mcout << "f.di=" << f.di << '\n';
50 mcout << "f.da=" << f.da << '\n';
51#endif
52#ifndef VISUAL_STUDIO
53 if (std::isnan(di) == 1) di = -DBL_MAX;
54 if (std::isnan(da) == 1) da = DBL_MAX;
55 if (std::isnan(f.di) == 1) f.di = -DBL_MAX;
56 if (std::isnan(f.da) == 1) f.da = DBL_MAX;
57#else
58 if (_isnan(di) == 1) di = -DBL_MAX;
59 if (_isnan(da) == 1) da = DBL_MAX;
60 if (_isnan(f.di) == 1) f.di = -DBL_MAX;
61 if (_isnan(f.da) == 1) f.da = DBL_MAX;
62#endif
63#ifdef POSSIBLE_FAILURE_ISNAN
64 if (!(di < d) && !(di >= d)) di = -DBL_MAX;
65 if (!(da < d) && !(da >= d)) da = DBL_MAX;
66 if (!(f.di < f.d) && !(f.di >= f.d)) f.di = -DBL_MAX;
67 if (!(f.da < f.d) && !(f.da >= f.d)) f.da = DBL_MAX;
68#endif
69#ifdef DEBUG_OPERATOR_MULT
70 DoubleAc temp(*this); // for debug
71#endif
72 d *= f.d;
73 if (di >= 0) {
74 if (f.di >= 0) {
75 di *= f.di;
76 da *= f.da;
77#ifdef DEBUG_OPERATOR_MULT_PRINT
78 mcout << "case 1\n";
79#endif
80 } else if (f.da >= 0) // assumed that f.di < 0
81 {
82 di = f.di * da;
83 da *= f.da;
84#ifdef DEBUG_OPERATOR_MULT_PRINT
85 mcout << "case 2\n";
86#endif
87 } else // assumed that f.da < 0
88 {
89 double ti = da * f.di;
90 da = di * f.da;
91 di = ti;
92#ifdef DEBUG_OPERATOR_MULT_PRINT
93 mcout << "case 3\n";
94#endif
95 }
96 } else if (da >= 0) // assumed that di < 0
97 {
98 if (f.di >= 0) {
99 di *= f.da;
100 da *= f.da;
101#ifdef DEBUG_OPERATOR_MULT_PRINT
102 mcout << "case 4\n";
103#endif
104 } else if (f.da >= 0) {
105 double ti = std::min(di * f.da, da * f.di);
106 da = std::max(di * f.di, da * f.da);
107 di = ti;
108#ifdef DEBUG_OPERATOR_MULT_PRINT
109 mcout << "case 5\n";
110#endif
111 } else {
112 double ti = da * f.di;
113 da = di * f.di;
114 di = ti;
115#ifdef DEBUG_OPERATOR_MULT_PRINT
116 mcout << "case 6\n";
117#endif
118 }
119 } else // assumed that di < 0 and da < 0
120 {
121 if (f.di >= 0) {
122 di *= f.da;
123 da *= f.di;
124#ifdef DEBUG_OPERATOR_MULT_PRINT
125 mcout << "case 7\n";
126#endif
127 } else if (f.da >= 0) {
128 double ti = di * f.da;
129 da = di * f.di;
130 di = ti;
131#ifdef DEBUG_OPERATOR_MULT_PRINT
132 mcout << "case 8\n";
133#endif
134 } else {
135 double ti = da * f.da;
136 da = di * f.di;
137 di = ti;
138#ifdef DEBUG_OPERATOR_MULT_PRINT
139 mcout << "case 8\n";
140#endif
141 }
142 }
143// mcout<<"d="<<d<<'\n';
144// mcout<<"di="<<di<<'\n';
145// mcout<<"da="<<da<<'\n';
146#ifdef CHECK_CORRECTNESS_AT_MULT
147 if (d < di) {
148 if (d - di == 0.0) // strange but at Sc. Linux 4.0 with -O2 this happens
149 goto mend; // but this precation does not cure!
150 funnw.ehdr(mcerr);
151 mcerr << "d < di at the end of computations\n";
152 mcerr << "This number:\n";
153 print(mcerr, 6);
154 mcerr << "Argument:\n";
155 f.print(mcerr, 6);
156 if (d > di) mcerr << "if(d > di) is also positive\n";
157 if (d == di) mcerr << "if(d == di) is also positive\n";
158 Iprintn(mcerr, d - di);
159#ifdef DEBUG_OPERATOR_MULT
160 // for debug:
161 mcerr << "old value:\n";
162 temp.print(mcerr, 6);
163#endif
164 spexit(mcerr);
165 }
166 if (d > da) {
167 if (d == da) // strange but at Sc. Linux 4.0 with O2 this happens
168 goto mend;
169 funnw.ehdr(mcerr);
170 mcerr << "d > di at the end of computations\n";
171 mcerr << "This number:\n";
172 print(mcerr, 6);
173 mcerr << "Argument:\n";
174 f.print(mcerr, 6);
175 if (d < da) mcerr << "if(d < da) is also positive\n";
176 if (d == da) mcerr << "if(d == da) is also positive\n";
177 Iprintn(mcerr, d - da);
178#ifdef DEBUG_OPERATOR_MULT
179 // for debug:
180 mcerr << "old value:\n";
181 temp.print(mcerr, 6);
182#endif
183 spexit(mcerr);
184 }
185mend:
186#endif
187 return *this;
188}
189
191 mfunnamep("DoubleAc& DoubleAc::operator/=(const DoubleAc& f)");
192 check_econd11(f.d, == 0, mcerr);
193 check_econd11(f.d, == 0.0, mcerr);
194#ifndef VISUAL_STUDIO
195 if (std::isnan(di) == 1) di = -DBL_MAX;
196 if (std::isnan(da) == 1) da = DBL_MAX;
197 if (std::isnan(f.di) == 1) f.di = -DBL_MAX;
198 if (std::isnan(f.da) == 1) f.da = DBL_MAX;
199#else
200 if (_isnan(di) == 1) di = -DBL_MAX;
201 if (_isnan(da) == 1) da = DBL_MAX;
202 if (_isnan(f.di) == 1) f.di = -DBL_MAX;
203 if (_isnan(f.da) == 1) f.da = DBL_MAX;
204#endif
205#ifdef POSSIBLE_FAILURE_ISNAN
206 if (!(di < d) && !(di >= d)) di = -DBL_MAX;
207 if (!(da < d) && !(da >= d)) da = DBL_MAX;
208 if (!(f.di < f.d) && !(f.di >= f.d)) f.di = -DBL_MAX;
209 if (!(f.da < f.d) && !(f.da >= f.d)) f.da = DBL_MAX;
210#endif
211
212 d /= f.d;
213 if (f.di < 0 && f.da > 0) {
214 di = -DBL_MAX;
215 da = DBL_MAX;
216 } else if (di >= 0) {
217 if (f.di > 0) {
218 di /= f.da;
219 da /= f.di;
220 } else if (f.di == 0) {
221 di /= f.da;
222 da = DBL_MAX;
223 } else {
224 if (f.da == 0) {
225 da = di / f.di;
226 di = -DBL_MAX;
227 } else {
228 double ti = da / f.da;
229 // mcout<<"d="<<d<<" ti="<<ti<<'\n';
230 da = di / f.di;
231 di = ti;
232 }
233 }
234 } else if (da >= 0) {
235 if (f.di > 0) {
236 di /= f.di;
237 da /= f.di;
238 } else if (f.di == 0) {
239 di = -DBL_MAX;
240 da = DBL_MAX;
241 } else {
242 if (f.da == 0) {
243 da = DBL_MAX;
244 di = -DBL_MAX;
245 } else {
246 double ti = da / f.da;
247 da = di / f.da;
248 di = ti;
249 }
250 }
251 } else {
252 // assumed da < 0
253 if (f.di > 0) {
254 di /= f.di;
255 da /= f.da;
256 } else if (f.di == 0) {
257 di = -DBL_MAX;
258 // check_econd11(f.da , == 0 , mcerr); // otherwise f.d == 0 which
259 // was already rejected
260 if (f.da == 0) {
261 funnw.ehdr(mcerr);
262 mcerr << "f.da == 0\n";
263 mcerr << "This means that f.d == 0 which should been already "
264 << "rejected.\n";
265 mcerr << "If the program reaches this point, this means that\n"
266 << "f.d is not between f.di and f.da, which is prohibited\n";
267 mcerr << "f :\n";
268 f.print(mcerr, 6);
269 spexit(mcerr);
270 }
271 da /= f.da;
272 } else {
273 if (f.da == 0) {
274 di = da / f.di;
275 da = DBL_MAX;
276 } else {
277 double ti = da / f.di;
278 da = di / f.da;
279 di = ti;
280 }
281 }
282 }
283#ifdef CHECK_CORRECTNESS_AT_MULT
284 if (d < di) {
285 funnw.ehdr(mcerr);
286 mcerr << "d < di at the end of computations\n";
287 mcerr << "This number:\n";
288 print(mcerr, 6);
289 mcerr << "Argument:\n";
290 f.print(mcerr, 6);
291 // for debug:
292 // mcerr<<"old value:\n";
293 // temp.print(mcerr, 6);
294 spexit(mcerr);
295 }
296 if (d > da) {
297 funnw.ehdr(mcerr);
298 mcerr << "d > di at the end of computations\n";
299 mcerr << "This number:\n";
300 print(mcerr, 6);
301 mcerr << "Argument:\n";
302 f.print(mcerr, 6);
303 // for debug:
304 // mcerr<<"old value:\n";
305 // temp.print(mcerr, 6);
306 spexit(mcerr);
307 }
308#endif
309 // check_econd12(d , < , di , mcerr);
310 // check_econd12(d , > , da , mcerr);
311 return *this;
312}
313
315 if (f.get() < 0) {
316 mcerr << "error in DoubleAc sqrt(const DoubleAc& f): f.get() < 0, f.get()="
317 << f.get() << '\n';
318 spexit(mcerr);
319 }
320 return DoubleAc(std::sqrt(double(f.get())),
321 std::sqrt(double(std::max(f.get_left_limit(), 0.0))),
322 std::sqrt(double(f.get_right_limit())));
323}
324
326 if (f.left_limit() >= 0.0) {
327 return DoubleAc(f.get() * f.get(), f.left_limit() * f.left_limit(),
328 f.right_limit() * f.right_limit());
329 } else if (f.right_limit() >= 0.0) {
330 double t = std::max(-f.left_limit(), f.right_limit());
331 return DoubleAc(f.get() * f.get(), 0.0, t * t);
332 }
333 return DoubleAc(f.get() * f.get(), f.right_limit() * f.right_limit(),
334 f.left_limit() * f.left_limit());
335}
336
337DoubleAc pow(const DoubleAc& f, double p) {
338 if (p == 1) return f;
339 if (p == 0) return DoubleAc(1.0);
340 if (p > 0) {
341 double d = std::pow(f.get(), p);
342 double di = std::pow(f.left_limit(), p);
343 double da = std::pow(f.right_limit(), p);
344 if (f.left_limit() >= 0.0) {
345 return DoubleAc(d, di, da);
346 } else if (f.right_limit() >= 0.0) {
347 if (di < 0.0)
348 return DoubleAc(d, di, da);
349 else // the power is even
350 return DoubleAc(d, 0.0, std::max(di, da));
351 } else {
352 if (di < 0.0)
353 return DoubleAc(d, di, da);
354 else // the power is even
355 return DoubleAc(d, da, di);
356 }
357 } else {
358 double d = std::pow(f.get(), -p);
359 double di = std::pow(f.left_limit(), -p);
360 double da = std::pow(f.right_limit(), -p);
361 if (f.left_limit() >= 0.0) {
362 return 1.0 / DoubleAc(d, di, da);
363 } else if (f.right_limit() >= 0.0) {
364 if (di < 0.0)
365 return 1.0 / DoubleAc(d, di, da);
366 else // the power is even
367 return 1.0 / DoubleAc(d, 0.0, std::max(di, da));
368 } else {
369 if (di < 0.0)
370 return 1.0 / DoubleAc(d, di, da);
371 else // the power is even
372 return 1.0 / DoubleAc(d, da, di);
373 }
374 }
375}
376
378 double d = std::exp(f.get());
379 double di = std::exp(f.left_limit());
380 double da = std::exp(f.right_limit());
381 return DoubleAc(d, di, da);
382}
383
385 // mcout<<"DoubleAc sin is starting, f="; f.print(mcout, 3);
386 double d = std::sin(f.get());
387 double di = std::sin(f.left_limit());
388 double da = std::sin(f.right_limit());
389 // Iprintn(mcout, d);
390 // Iprintn(mcout, di);
391 // Iprintn(mcout, da);
392 long n = left_round(f.get() / M_PI + 0.5);
393 long ni = left_round(f.left_limit() / M_PI + 0.5);
394 long na = left_round(f.right_limit() / M_PI + 0.5);
395 // Iprintn(mcout, n);
396 // Iprintn(mcout, ni);
397 // Iprintn(mcout, na);
398 if (n % 2 == 0) {
399 // Even number.
400 if (ni < n) {
401 di = -1.0;
402 da = std::max(di, da);
403 if (na > n) {
404 da = 1.0;
405 }
406 } else if (na > n) {
407 da = 1.0;
408 di = std::min(di, da);
409 }
410 } else {
411 // Odd number.
412 double temp = di;
413 di = da;
414 da = temp;
415 if (ni < n) {
416 da = 1.0;
417 di = std::min(di, da);
418 if (na > n) {
419 di = -1.0;
420 }
421 } else if (na > n) {
422 di = -1.0;
423 da = std::max(di, da);
424 }
425 }
426 // Iprintn(mcout, d);
427 // Iprintn(mcout, di);
428 // Iprintn(mcout, da);
429 return DoubleAc(d, di, da);
430}
431
433 double d = std::cos(f.get());
434 double di = std::cos(f.left_limit());
435 double da = std::cos(f.right_limit());
436 long n = left_round(f.get() / M_PI - 1.0);
437 long ni = left_round(f.left_limit() / M_PI - 1.0);
438 long na = left_round(f.right_limit() / M_PI - 1.0);
439 if (n % 2 == 0) {
440 // Even number.
441 if (ni < n) {
442 di = -1.0;
443 da = std::max(di, da);
444 if (na > n) {
445 da = 1.0;
446 }
447 } else if (na > n) {
448 da = 1.0;
449 di = std::min(di, da);
450 }
451 } else {
452 // Odd number.
453 double temp = di;
454 di = da;
455 da = temp;
456 if (ni < n) {
457 da = 1.0;
458 di = std::min(di, da);
459 if (na > n) {
460 di = -1.0;
461 }
462 } else if (na > n) {
463 di = -1.0;
464 da = std::max(di, da);
465 }
466 }
467 return DoubleAc(d, di, da);
468}
469
471 if (fabs(f.get()) > 1) {
472 mcerr << "ERROR in inline DoubleAc asin(const DoubleAc& f):\n"
473 << "fabs(f.get()) > 1: f.get()=" << f.get() << '\n';
474 spexit(mcerr);
475 }
476 double d = std::asin(f.get());
477 double di;
478 if (f.left_limit() < -1.0)
479 di = std::asin(-1.0);
480 else
481 di = std::asin(f.left_limit());
482 double da;
483 if (f.right_limit() > 1.0)
484 da = std::asin(1.0);
485 else
486 da = std::asin(f.right_limit());
487 return DoubleAc(d, di, da);
488}
489
491 if (fabs(f.get()) > 1) {
492 mcerr << "ERROR in inline DoubleAc acos(const DoubleAc& f):\n"
493 << "fabs(f.get()) > 1: f.get()=" << f.get() << '\n';
494 spexit(mcerr);
495 }
496 double d = std::acos(f.get());
497 double da;
498 if (f.left_limit() < -1.0)
499 da = std::acos(-1.0);
500 else
501 da = std::acos(f.left_limit());
502 double di;
503 if (f.right_limit() > 1.0)
504 di = std::acos(1.0);
505 else
506 di = std::acos(f.right_limit());
507 return DoubleAc(d, di, da);
508}
509
510void DoubleAc::print(std::ostream& file, int l) const {
511 if (l <= 0) return;
512 if (l == 1) {
513 file << d;
514 } else if (l == 2) {
515 file << d << " [ " << di << " , " << da << " ] ";
516 } else if (l == 3) {
517 file << d;
518 int t = file.precision(2);
519 file << " [" << std::setw(8) << d - di << "," << std::setw(8) << da - d
520 << "] ";
521 file.precision(t);
522 } else if (l == 4) {
523 file << d << " [ " << di << " , " << da << " ] \n";
524 } else if (l == 5) {
525 file << d;
526 int t = file.precision(2);
527 file << " [" << std::setw(8) << d - di << "," << std::setw(8) << da - d
528 << "] \n";
529 file.precision(t);
530 } else {
531 int t = file.precision(16);
532 file << "DoubleAc: d=" << std::setw(20) << d << " di=" << std::setw(20)
533 << di << " da=" << std::setw(20) << da << '\n';
534 file.precision(t);
535 }
536}
537
538DoubleAc pow(const DoubleAc& /*f*/, const DoubleAc& /*p*/) {
539 mcerr << "ERROR in inline DoubleAc pow(const DoubleAc& f, const DoubleAc& "
540 "p):\n";
541 mcerr << "not implemented yet\n";
542 spexit(mcerr);
543 return 0.0;
544}
545
546std::ostream& operator<<(std::ostream& file, const DoubleAc& f) {
547 f.print(file, 1);
548 return file;
549}
550
551}
#define check_econd11(a, signb, stream)
Definition: FunNameStack.h:155
#define check_econd12a(a, sign, b, add, stream)
Definition: FunNameStack.h:181
#define mfunnamep(string)
Definition: FunNameStack.h:49
#define spexit(stream)
Definition: FunNameStack.h:256
#define mfunname(string)
Definition: FunNameStack.h:45
double get_right_limit(void) const
Definition: DoubleAc.h:82
DoubleAc & operator/=(DoubleAc f)
Definition: DoubleAc.cpp:190
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
double get_left_limit(void) const
Definition: DoubleAc.h:78
double left_limit(void) const
Definition: DoubleAc.h:79
double right_limit(void) const
Definition: DoubleAc.h:83
DoubleAc & operator*=(DoubleAc f)
Definition: DoubleAc.cpp:41
Definition: BGMesh.cpp:6
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
std::ostream & operator<<(std::ostream &file, const BGMesh &bgm)
Definition: BGMesh.cpp:37
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:337
long left_round(double f)
Definition: minmax.h:45
DoubleAc exp(const DoubleAc &f)
Definition: DoubleAc.cpp:377
DoubleAc square(const DoubleAc &f)
Definition: DoubleAc.cpp:325
DoubleAc acos(const DoubleAc &f)
Definition: DoubleAc.cpp:490
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384
DoubleAc asin(const DoubleAc &f)
Definition: DoubleAc.cpp:470
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314
#define mcout
Definition: prstream.h:126
#define mcerr
Definition: prstream.h:128
#define Iprintn(file, name)
Definition: prstream.h:205