1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
//
// A rust binding for the GSL library by Guillaume Gomez (guillaume1.gomez@gmail.com)
//

#[derive(Clone, PartialEq, PartialOrd, Debug, Copy)]
#[repr(C)]
pub enum Mode {
    PrecDouble,
    PrecSingle,
    PrecApprox,
}

/// A type for results generated by GSL functions where `Err` is `enums::Value`.
pub type GSLResult<T> = ::std::result::Result<T, Value>;

impl ::std::convert::From<Value> for GSLResult<()> {
    fn from(v: Value) -> Self {
        match v {
            Value::Success => Ok(()),
            e => { Err(e) },
        }
    }
}

#[derive(Clone, PartialEq, PartialOrd, Debug, Copy)]
#[repr(C)]
pub enum Value {
    Success = 0,
    Failure = -1,
    /// iteration has not converged
    Continue = -2,
    /// input domain error, e.g sqrt(-1)
    Domain = 1,
    /// output range error, e.g. exp(1e100)
    Range = 2,
    /// invalid pointer
    Fault = 3,
    /// invalid argument supplied by user
    Invalid = 4,
    /// generic failure
    Failed = 5,
    /// factorization failed
    Factorization = 6,
    /// sanity check failed - shouldn't happen
    Sanity = 7,
    /// malloc failed
    NoMemory = 8,
    /// problem with user-supplied function
    BadFunction = 9,
    /// iterative process is out of control
    RunAway = 10,
    /// exceeded max number of iterations
    MaxIteration = 11,
    /// tried to divide by zero
    ZeroDiv = 12,
    /// user specified an invalid tolerance
    BadTolerance = 13,
    /// failed to reach the specified tolerance
    Tolerance = 14,
    /// underflow
    UnderFlow = 15,
    /// overflow
    OverFlow = 16,
    /// loss of accuracy
    Loss = 17,
    /// failed because of roundoff error
    Round = 18,
    /// matrix, vector lengths are not conformant
    BadLength = 19,
    /// matrix not square
    NotSquare = 20,
    /// apparent singularity detected
    Singularity = 21,
    /// integral or series is divergent
    Diverge = 22,
    /// requested feature is not supported by the hardware
    Unsupported = 23,
    /// requested feature not (yet) implemented
    Unimplemented = 24,
    /// cache limit exceeded
    Cache = 25,
    /// table limit exceeded
    Table = 26,
    /// iteration is not making progress towards solution
    NoProgress = 27,
    /// jacobian evaluations are not improving the solution
    NoProgressJacobian = 28,
    /// cannot reach the specified tolerance in F
    ToleranceF = 29,
    /// cannot reach the specified tolerance in X
    ToleranceX = 30,
    /// cannot reach the specified tolerance in gradient
    ToleranceG = 31,
    /// cannot reach the specified tolerance in gradient
    EOF = 32,
}

impl Value {
    pub fn wrap(v: i32) -> Value {
        match v {
            -2 => Value::Continue,
            -1 => Value::Failure,
            0  => Value::Success,
            1  => Value::Domain,
            2  => Value::Range,
            3  => Value::Fault,
            4  => Value::Invalid,
            5  => Value::Failed,
            6  => Value::Factorization,
            7  => Value::Sanity,
            8  => Value::NoMemory,
            9  => Value::BadFunction,
            10 => Value::RunAway,
            11 => Value::MaxIteration,
            12 => Value::ZeroDiv,
            13 => Value::BadTolerance,
            14 => Value::Tolerance,
            15 => Value::UnderFlow,
            16 => Value::OverFlow,
            17 => Value::Loss,
            18 => Value::Round,
            19 => Value::BadLength,
            20 => Value::NotSquare,
            21 => Value::Singularity,
            22 => Value::Diverge,
            23 => Value::Unsupported,
            24 => Value::Unimplemented,
            25 => Value::Cache,
            26 => Value::Table,
            27 => Value::NoProgress,
            28 => Value::NoProgressJacobian,
            29 => Value::ToleranceF,
            30 => Value::ToleranceX,
            31 => Value::ToleranceG,
            32 => Value::EOF,
            x  => panic!("Unknow value for Value enum: '{}'", x),
        }
    }
}

#[derive(Clone, PartialEq, PartialOrd, Debug, Copy)]
#[repr(C)]
pub enum EigenSort {
    /// ascending order in numerical value
    ValAsc,
    /// descending order in numerical value
    ValDesc,
    /// ascending order in magnitude
    AbsAsc,
    /// descending order in magnitude
    AbsDesc,
}

#[derive(Clone, PartialEq, PartialOrd, Debug, Copy)]
#[repr(C)]
/// this gives the sign in the formula
///
/// h(f) = \sum x(t) exp(+/- 2 pi i f t)
///
/// where - is the forward transform direction and + the inverse direction
pub enum FftDirection {
    Forward = -1,
    Backward = 1,
}

#[derive(Clone, PartialEq, PartialOrd, Debug, Copy)]
#[repr(C)]
/// The low-level integration rules in QUADPACK are identified by small integers (1-6). We'll use symbolic constants to refer to them.
pub enum GaussKonrodRule {
    /// 15 point Gauss-Kronrod rule
    Gauss15 = 1,
    /// 21 point Gauss-Kronrod rule
    Gauss21 = 2,
    /// 31 point Gauss-Kronrod rule
    Gauss31 = 3,
    /// 41 point Gauss-Kronrod rule
    Gauss41 = 4,
    /// 51 point Gauss-Kronrod rule
    Gauss51 = 5,
    /// 61 point Gauss-Kronrod rule
    Gauss61 = 6,
}

#[derive(Clone, PartialEq, PartialOrd, Debug, Copy)]
#[repr(C)]
/// Used by workspace for QAWO integrator
pub enum IntegrationQawo {
    Cosine,
    Sine,
}

#[derive(Clone, PartialEq, PartialOrd, Debug, Copy)]
#[repr(C)]
/// Used by VegasMonteCarlo struct
///
/// The possible choices are GSL_VEGAS_MODE_IMPORTANCE, GSL_VEGAS_MODE_
/// STRATIFIED, GSL_VEGAS_MODE_IMPORTANCE_ONLY. This determines whether vegas
/// will use importance sampling or stratified sampling, or whether it can pick on
/// its own. In low dimensions vegas uses strict stratified sampling (more precisely,
/// stratified sampling is chosen if there are fewer than 2 bins per box).
pub enum VegasMode {
    Importance = 1,
    ImportanceOnly = 0,
    Stratified = -1,
}

#[derive(Clone, PartialEq, PartialOrd, Debug, Copy)]
#[repr(C)]
/// Possible return values for an hadjust() evolution method for ordinary differential equations
pub enum ODEiv {
    /// step was increased
    Inc = 1,
    /// step unchanged
    Nil = 0,
    /// step decreased
    Dec = -1,
}

#[derive(Clone, PartialEq, PartialOrd, Debug, Copy)]
#[repr(C)]
pub enum WaveletDirection {
    Forward = 1,
    Backward = -1,
}