Garfield++ v2r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Heed::DoubleAc Class Reference

#include <DoubleAc.h>

Public Member Functions

 DoubleAc (void)
 
 DoubleAc (const DoubleAc &f)
 
 DoubleAc (double f)
 
 DoubleAc (float f)
 
 DoubleAc (long f)
 
 DoubleAc (int f)
 
 DoubleAc (double f, double ffmin, double ffmax)
 
 DoubleAc (double f, double relative_prec)
 
DoubleAcoperator= (const DoubleAc &f)
 
DoubleAcoperator= (double f)
 
DoubleAcoperator= (float f)
 
DoubleAcoperator= (long f)
 
DoubleAcoperator= (int f)
 
 operator double (void) const
 
double get (void) const
 
double get_low_limit (void) const
 
double get_left_limit (void) const
 
double left_limit (void) const
 
double get_min_limit (void) const
 
double get_high_limit (void) const
 
double get_right_limit (void) const
 
double right_limit (void) const
 
double get_max_limit (void) const
 
double get_accuracy (void) const
 
DoubleAcoperator+= (const DoubleAc &f)
 
DoubleAcoperator+= (double f)
 
DoubleAcoperator+= (float f)
 
DoubleAcoperator+= (long f)
 
DoubleAcoperator+= (int f)
 
DoubleAcoperator-= (const DoubleAc &f)
 
DoubleAcoperator-= (double f)
 
DoubleAcoperator-= (float f)
 
DoubleAcoperator-= (long f)
 
DoubleAcoperator-= (int f)
 
DoubleAcoperator*= (DoubleAc f)
 
DoubleAcoperator*= (double f)
 
DoubleAcoperator*= (float f)
 
DoubleAcoperator*= (long f)
 
DoubleAcoperator*= (int f)
 
DoubleAcoperator/= (DoubleAc f)
 
DoubleAcoperator/= (double f)
 
DoubleAcoperator/= (float f)
 
DoubleAcoperator/= (long f)
 
DoubleAcoperator/= (int f)
 
void print (std::ostream &file, int l=1) const
 

Friends

void change_sign (DoubleAc &f)
 

Detailed Description

"Double with accuracy". Really this is original implementation of interval computations. This is the way of watching for the numerical errors by establishing lower and upper limit of each numerical value. The central the most probable value is also watched for.

Definition at line 44 of file DoubleAc.h.

Constructor & Destructor Documentation

◆ DoubleAc() [1/8]

Heed::DoubleAc::DoubleAc ( void  )
inline

Definition at line 49 of file DoubleAc.h.

49 {
50 di = 0.0;
51 d = 0.0;
52 da = 0.0;
53 }

Referenced by operator+=(), and operator-=().

◆ DoubleAc() [2/8]

Heed::DoubleAc::DoubleAc ( const DoubleAc f)
inline

Definition at line 54 of file DoubleAc.h.

54 {
55 di = f.di;
56 d = f.d;
57 da = f.da;
58 }

◆ DoubleAc() [3/8]

Heed::DoubleAc::DoubleAc ( double  f)
inline

Definition at line 130 of file DoubleAc.h.

130 {
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}
const double one_plus_def_dbl_prec
Definition: DoubleAc.h:31

◆ DoubleAc() [4/8]

Heed::DoubleAc::DoubleAc ( float  f)
inline

Definition at line 189 of file DoubleAc.h.

189 {
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}
const double one_plus_def_flt_prec
Definition: DoubleAc.h:35
const double one_minus_def_flt_prec
Definition: DoubleAc.h:36

◆ DoubleAc() [5/8]

Heed::DoubleAc::DoubleAc ( long  f)
inline

Definition at line 218 of file DoubleAc.h.

218 {
219 d = f;
220 di = f;
221 da = f;
222}

◆ DoubleAc() [6/8]

Heed::DoubleAc::DoubleAc ( int  f)
inline

Definition at line 231 of file DoubleAc.h.

231 {
232 d = f;
233 di = f;
234 da = f;
235}

◆ DoubleAc() [7/8]

Heed::DoubleAc::DoubleAc ( double  f,
double  ffmin,
double  ffmax 
)

Definition at line 18 of file DoubleAc.cpp.

18 {
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}
#define check_econd12a(a, sign, b, add, stream)
Definition: FunNameStack.h:181
#define mfunname(string)
Definition: FunNameStack.h:45
#define mcerr
Definition: prstream.h:128

◆ DoubleAc() [8/8]

Heed::DoubleAc::DoubleAc ( double  f,
double  relative_prec 
)

Definition at line 27 of file DoubleAc.cpp.

27 {
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}
#define check_econd11(a, signb, stream)
Definition: FunNameStack.h:155

Member Function Documentation

◆ get()

double Heed::DoubleAc::get ( void  ) const
inline

◆ get_accuracy()

double Heed::DoubleAc::get_accuracy ( void  ) const
inline

Definition at line 85 of file DoubleAc.h.

85{ return 0.5 * (da - di); }

◆ get_high_limit()

double Heed::DoubleAc::get_high_limit ( void  ) const
inline

Definition at line 81 of file DoubleAc.h.

81{ return da; }

◆ get_left_limit()

double Heed::DoubleAc::get_left_limit ( void  ) const
inline

Definition at line 78 of file DoubleAc.h.

78{ return di; }

Referenced by Heed::operator-(), and Heed::sqrt().

◆ get_low_limit()

double Heed::DoubleAc::get_low_limit ( void  ) const
inline

Definition at line 77 of file DoubleAc.h.

77{ return di; }

◆ get_max_limit()

double Heed::DoubleAc::get_max_limit ( void  ) const
inline

Definition at line 84 of file DoubleAc.h.

84{ return da; }

◆ get_min_limit()

double Heed::DoubleAc::get_min_limit ( void  ) const
inline

Definition at line 80 of file DoubleAc.h.

80{ return di; }

◆ get_right_limit()

double Heed::DoubleAc::get_right_limit ( void  ) const
inline

Definition at line 82 of file DoubleAc.h.

82{ return da; }

Referenced by Heed::operator-(), and Heed::sqrt().

◆ left_limit()

double Heed::DoubleAc::left_limit ( void  ) const
inline

◆ operator double()

Heed::DoubleAc::operator double ( void  ) const
inline

Definition at line 75 of file DoubleAc.h.

75{ return d; }

◆ operator*=() [1/5]

DoubleAc & Heed::DoubleAc::operator*= ( double  f)
inline

Definition at line 301 of file DoubleAc.h.

301 {
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}

◆ operator*=() [2/5]

DoubleAc & Heed::DoubleAc::operator*= ( DoubleAc  f)

Definition at line 41 of file DoubleAc.cpp.

41 {
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}
#define mfunnamep(string)
Definition: FunNameStack.h:49
#define spexit(stream)
Definition: FunNameStack.h:256
DoubleAc(void)
Definition: DoubleAc.h:49
void print(std::ostream &file, int l=1) const
Definition: DoubleAc.cpp:510
#define mcout
Definition: prstream.h:126
#define Iprintn(file, name)
Definition: prstream.h:205

◆ operator*=() [3/5]

DoubleAc & Heed::DoubleAc::operator*= ( float  f)
inline

Definition at line 313 of file DoubleAc.h.

313 {
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}

◆ operator*=() [4/5]

DoubleAc & Heed::DoubleAc::operator*= ( int  f)
inline

Definition at line 337 of file DoubleAc.h.

337 {
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}

◆ operator*=() [5/5]

DoubleAc & Heed::DoubleAc::operator*= ( long  f)
inline

Definition at line 325 of file DoubleAc.h.

325 {
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}

◆ operator+=() [1/5]

DoubleAc & Heed::DoubleAc::operator+= ( const DoubleAc f)
inline

Definition at line 243 of file DoubleAc.h.

243 {
244 di += f.di;
245 d += f.d;
246 da += f.da;
247 return *this;
248}

◆ operator+=() [2/5]

DoubleAc & Heed::DoubleAc::operator+= ( double  f)
inline

Definition at line 249 of file DoubleAc.h.

249 {
250 *this += DoubleAc(f);
251 return *this;
252}

◆ operator+=() [3/5]

DoubleAc & Heed::DoubleAc::operator+= ( float  f)
inline

Definition at line 254 of file DoubleAc.h.

254 {
255 *this += DoubleAc(f);
256 return *this;
257}

◆ operator+=() [4/5]

DoubleAc & Heed::DoubleAc::operator+= ( int  f)
inline

Definition at line 265 of file DoubleAc.h.

265 {
266 di += f;
267 d += f;
268 da += f;
269 return *this;
270}

◆ operator+=() [5/5]

DoubleAc & Heed::DoubleAc::operator+= ( long  f)
inline

Definition at line 259 of file DoubleAc.h.

259 {
260 di += f;
261 d += f;
262 da += f;
263 return *this;
264}

◆ operator-=() [1/5]

DoubleAc & Heed::DoubleAc::operator-= ( const DoubleAc f)
inline

Definition at line 272 of file DoubleAc.h.

272 {
273 di -= f.da;
274 d -= f.d;
275 da -= f.di;
276 return *this;
277}

◆ operator-=() [2/5]

DoubleAc & Heed::DoubleAc::operator-= ( double  f)
inline

Definition at line 278 of file DoubleAc.h.

278 {
279 *this -= DoubleAc(f);
280 return *this;
281}

◆ operator-=() [3/5]

DoubleAc & Heed::DoubleAc::operator-= ( float  f)
inline

Definition at line 283 of file DoubleAc.h.

283 {
284 *this -= DoubleAc(f);
285 return *this;
286}

◆ operator-=() [4/5]

DoubleAc & Heed::DoubleAc::operator-= ( int  f)
inline

Definition at line 294 of file DoubleAc.h.

294 {
295 di -= f;
296 d -= f;
297 da -= f;
298 return *this;
299}

◆ operator-=() [5/5]

DoubleAc & Heed::DoubleAc::operator-= ( long  f)
inline

Definition at line 288 of file DoubleAc.h.

288 {
289 di -= f;
290 d -= f;
291 da -= f;
292 return *this;
293}

◆ operator/=() [1/5]

DoubleAc & Heed::DoubleAc::operator/= ( double  f)
inline

Definition at line 350 of file DoubleAc.h.

350 {
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}

◆ operator/=() [2/5]

DoubleAc & Heed::DoubleAc::operator/= ( DoubleAc  f)

Definition at line 190 of file DoubleAc.cpp.

190 {
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}

◆ operator/=() [3/5]

DoubleAc & Heed::DoubleAc::operator/= ( float  f)
inline

Definition at line 367 of file DoubleAc.h.

367 {
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}

◆ operator/=() [4/5]

DoubleAc & Heed::DoubleAc::operator/= ( int  f)
inline

Definition at line 401 of file DoubleAc.h.

401 {
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}

◆ operator/=() [5/5]

DoubleAc & Heed::DoubleAc::operator/= ( long  f)
inline

Definition at line 384 of file DoubleAc.h.

384 {
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}

◆ operator=() [1/5]

DoubleAc & Heed::DoubleAc::operator= ( const DoubleAc f)
inline

Definition at line 65 of file DoubleAc.h.

65 {
66 d = f.d;
67 di = f.di;
68 da = f.da;
69 return *this;
70 }

◆ operator=() [2/5]

DoubleAc & Heed::DoubleAc::operator= ( double  f)
inline

Definition at line 160 of file DoubleAc.h.

160 {
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}

◆ operator=() [3/5]

DoubleAc & Heed::DoubleAc::operator= ( float  f)
inline

Definition at line 203 of file DoubleAc.h.

203 {
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}

◆ operator=() [4/5]

DoubleAc & Heed::DoubleAc::operator= ( int  f)
inline

Definition at line 236 of file DoubleAc.h.

236 {
237 d = f;
238 di = f;
239 da = f;
240 return *this;
241}

◆ operator=() [5/5]

DoubleAc & Heed::DoubleAc::operator= ( long  f)
inline

Definition at line 224 of file DoubleAc.h.

224 {
225 d = f;
226 di = f;
227 da = f;
228 return *this;
229}

◆ print()

void Heed::DoubleAc::print ( std::ostream &  file,
int  l = 1 
) const

Definition at line 510 of file DoubleAc.cpp.

510 {
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}

Referenced by operator*=(), operator/=(), and Heed::operator<<().

◆ right_limit()

double Heed::DoubleAc::right_limit ( void  ) const
inline

Definition at line 83 of file DoubleAc.h.

83{ return da; }

Referenced by Heed::acos(), Heed::asin(), Heed::cos(), Heed::exp(), Heed::fabs(), Heed::pow(), Heed::sin(), and Heed::square().

Friends And Related Function Documentation

◆ change_sign

void change_sign ( DoubleAc f)
friend

Definition at line 424 of file DoubleAc.h.

424 {
425 f.d = -f.d;
426 double temp = f.di;
427 f.di = -f.da;
428 f.da = -temp;
429}

The documentation for this class was generated from the following files: