summaryrefslogtreecommitdiff
path: root/gsl-1.9/multimin/diff.c
blob: 218cc7fb0d3ffa9c79f28823582602352d4e1ba1 (plain)
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
/* multimin/diff.c
 * 
 * Copyright (C) 2000 David Morrison
 * 
 * This program is free software; you can redistribute it and/or modify
 * it under the terms of the GNU General Public License as published by
 * the Free Software Foundation; either version 2 of the License, or (at
 * your option) any later version.
 * 
 * This program is distributed in the hope that it will be useful, but
 * WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
 * General Public License for more details.
 * 
 * You should have received a copy of the GNU General Public License
 * along with this program; if not, write to the Free Software
 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
 */

#include <config.h>
#include <gsl/gsl_multimin.h>

int
gsl_multimin_diff (const gsl_multimin_function * f,
                   const gsl_vector * x, gsl_vector * g)
{
  size_t i, n = f->n;

  double h = GSL_SQRT_DBL_EPSILON;


  gsl_vector * x1 = gsl_vector_alloc (n);  /* FIXME: pass as argument */

  gsl_vector_memcpy (x1, x);

  for (i = 0; i < n; i++)
    {
      double fl, fh;
  
      double xi = gsl_vector_get (x, i);
      double dx = fabs(xi) * h;

      if (dx == 0.0) dx = h;

      gsl_vector_set (x1, i, xi + dx);
      fh = GSL_MULTIMIN_FN_EVAL(f, x1);

      gsl_vector_set (x1, i, xi - dx);
      fl = GSL_MULTIMIN_FN_EVAL(f, x1);

      gsl_vector_set (x1, i, xi);
      gsl_vector_set (g, i, (fh - fl) / (2.0 * dx));
    }

  gsl_vector_free (x1);

  return GSL_SUCCESS;
}