Dune Core Modules (unstable)

gmpfield.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5#ifndef DUNE_GMPFIELD_HH
6#define DUNE_GMPFIELD_HH
7
12#include <iostream>
13#include <string>
14#include <type_traits>
15
16#if HAVE_GMP || DOXYGEN
17
18#include <gmpxx.h>
19
20#include <dune/common/math.hh>
23
24namespace Dune
25{
26
31 template< unsigned int precision >
33 : public mpf_class
34 {
35 typedef mpf_class Base;
36
37 public:
40 : Base(0,precision)
41 {}
42
46 GMPField ( const char* str )
47 : Base(str,precision)
48 {}
49
53 GMPField ( const std::string& str )
54 : Base(str,precision)
55 {}
56
59 template< class T,
60 typename EnableIf = typename std::enable_if<
61 std::is_convertible<T, mpf_class>::value>::type
62 >
63 GMPField ( const T &v )
64 : Base( v,precision )
65 {}
66
67 // type conversion operators
68 operator double () const
69 {
70 return this->get_d();
71 }
72
73 };
74
75 template <unsigned int precision>
76 struct IsNumber<GMPField<precision>>
77 : public std::integral_constant<bool, true> {
78 };
79
80 template< unsigned int precision1, unsigned int precision2 >
81 struct PromotionTraits<GMPField<precision1>, GMPField<precision2>>
82 {
83 typedef GMPField<(precision1 > precision2 ? precision1 : precision2)> PromotedType;
84 };
85
86 template< unsigned int precision >
87 struct PromotionTraits<GMPField<precision>,GMPField<precision>>
88 {
89 typedef GMPField<precision> PromotedType;
90 };
91
92 template< unsigned int precision, class T >
93 struct PromotionTraits<GMPField<precision>, T>
94 {
95 typedef GMPField<precision> PromotedType;
96 };
97
98 template< class T, unsigned int precision >
99 struct PromotionTraits<T, GMPField<precision>>
100 {
101 typedef GMPField<precision> PromotedType;
102 };
103
104 template< unsigned int precision >
105 struct MathematicalConstants<GMPField<precision>>
106 {
107 static_assert(precision < 3319u, "Mathematical constants for GMPField defined up to a precision 3318."); // 1000 digits10
108
109 using T = GMPField<precision>;
110 static const T e ()
111 {
112 static const T e = T("2.7182818284590452353602874713526624977572470936999595749669676277240766303535475945713821785251664274274663919320030599218174135966290435729003342952605956307381323286279434907632338298807531952510190115738341879307021540891499348841675092447614606680822648001684774118537423454424371075390777449920695517027618386062613313845830007520449338265602976067371132007093287091274437470472306969772093101416928368190255151086574637721112523897844250569536967707854499699679468644549059879316368892300987931277361782154249992295763514822082698951936680331825288693984964651058209392398294887933203625094431173012381970684161403970198376793206832823764648042953118023287825098194558153017567173613320698112509961818815930416903515988885193458072738667385894228792284998920868058257492796104841984443634632449684875602336248270419786232090021609902353043699418491463140934317381436405462531520961836908887070167683964243781405927145635490613031072085103837505101157477041718986106873969655212671546889570350");
113 return e;
114 }
115
116 static const T pi ()
117 {
118 static const T pi = T("3.1415926535897932384626433832795028841971693993751058209749445923078164062862089986280348253421170679821480865132823066470938446095505822317253594081284811174502841027019385211055596446229489549303819644288109756659334461284756482337867831652712019091456485669234603486104543266482133936072602491412737245870066063155881748815209209628292540917153643678925903600113305305488204665213841469519415116094330572703657595919530921861173819326117931051185480744623799627495673518857527248912279381830119491298336733624406566430860213949463952247371907021798609437027705392171762931767523846748184676694051320005681271452635608277857713427577896091736371787214684409012249534301465495853710507922796892589235420199561121290219608640344181598136297747713099605187072113499999983729780499510597317328160963185950244594553469083026425223082533446850352619311881710100031378387528865875332083814206171776691473035982534904287554687311595628638823537875937519577818577805321712268066130019278766111959092164201");
119 return pi;
120 }
121 };
122}
123
124#endif // HAVE_GMP
125
126#endif // #ifndef DUNE_GMPFIELD_HH
Number class for high precision floating point number using the GMP library mpf_class implementation.
Definition: gmpfield.hh:34
GMPField(const T &v)
initialize from a compatible scalar type
Definition: gmpfield.hh:63
GMPField()
Definition: gmpfield.hh:39
GMPField(const std::string &str)
initialize from a string
Definition: gmpfield.hh:53
GMPField(const char *str)
initialize from a string
Definition: gmpfield.hh:46
Some useful basic math stuff.
Dune namespace.
Definition: alignedallocator.hh:13
Compute type of the result of an arithmetic operation involving two different number types.
static const Field e()
Euler's number.
Definition: math.hh:38
static const Field pi()
Archimedes' constant.
Definition: math.hh:48
Traits for type conversions and type information.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)