Octopus
oct_gsl_f.c
Go to the documentation of this file.
1/*
2 Copyright (C) 2002 M. Marques, A. Castro, A. Rubio, G. Bertsch, M. Oliveira
3
4 This program is free software; you can redistribute it and/or modify
5 it under the terms of the GNU General Public License as published by
6 the Free Software Foundation; either version 2, or (at your option)
7 any later version.
8
9 This program is distributed in the hope that it will be useful,
10 but WITHOUT ANY WARRANTY; without even the implied warranty of
11 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 GNU General Public License for more details.
13
14 You should have received a copy of the GNU General Public License
15 along with this program; if not, write to the Free Software
16 Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17 02110-1301, USA.
18
19*/
20
21#include <config.h>
22
23#include <math.h>
24#include <stdio.h>
25#include <stdlib.h>
26
27#include <gsl/gsl_combination.h>
28#include <gsl/gsl_randist.h>
29#include <gsl/gsl_rng.h>
30#include <gsl/gsl_sf.h>
31
32#include "string_f.h"
33
34#include <fortran_types.h>
35
36/* Numerical threshold for oct_bessel_k0 and oct_bessel_k1 */
37#define BESSEL_K_THRES 1.0e2
38
39/* ---------------------- Interface to GSL functions ------------------------ */
40
41/* Special Functions */
42double FC_FUNC_(oct_gamma, OCT_GAMMA)(const double *x) {
43 return gsl_sf_gamma(*x);
44}
45
46double FC_FUNC_(oct_incomplete_gamma, OCT_INCOMPLETE_GAMMA)(const double *a,
47 const double *x) {
48 return gsl_sf_gamma_inc_Q(*a, *x);
49}
50
51double FC_FUNC_(oct_sph_bessel, OCT_SPH_BESSEL)(const fint *l,
52 const double *x) {
53 if (*l == 0) {
54 return gsl_sf_bessel_j0(*x);
55 } else if (*l == 1) {
56 return gsl_sf_bessel_j1(*x);
57 } else if (*l == 2) {
58 return gsl_sf_bessel_j2(*x);
59 } else {
60 return gsl_sf_bessel_jl(*l, *x);
61 }
62}
63
64double FC_FUNC_(oct_bessel, OCT_BESSEL)(const fint *n, const double *x) {
65 return gsl_sf_bessel_Jn(*n, *x);
68double FC_FUNC_(oct_bessel_in, OCT_BESSEL_IN)(const fint *n, const double *x) {
69 return gsl_sf_bessel_In(*n, *x);
70}
72double FC_FUNC_(oct_bessel_j0, OCT_BESSEL_J0)(const double *x) {
73 return gsl_sf_bessel_J0(*x);
76double FC_FUNC_(oct_bessel_j1, OCT_BESSEL_J1)(const double *x) {
77 return gsl_sf_bessel_J1(*x);
80double FC_FUNC_(oct_bessel_k0, OCT_BESSEL_K0)(const double *x) {
81 if (*x > BESSEL_K_THRES) {
82 return 0.0e0;
83 } else {
84 return gsl_sf_bessel_K0(*x);
85 }
88double FC_FUNC_(oct_bessel_k1, OCT_BESSEL_K1)(const double *x) {
89 if (*x > BESSEL_K_THRES) {
90 return 0.0e0;
91 } else {
92 return gsl_sf_bessel_K1(*x);
93 }
94}
95
96/* the GSL headers specify double x, not const double x */
97double FC_FUNC_(oct_erfc, OCT_ERFC)(const double *x) {
98 /* avoid floating invalids in the asymptotic limit */
99 if (*x > 20.0)
100 return 0.0;
101 if (*x < -20.0)
102 return 2.0;
103 /* otherwise call gsl */
104 return gsl_sf_erfc(*x);
105}
107/* the GSL headers specify double x, not const double x */
108double FC_FUNC_(oct_erf, OCT_ERF)(const double *x) {
109 /* avoid floating invalids in the asymptotic limit */
110 if (*x > 20.0)
111 return 1.0;
112 if (*x < -20.0)
113 return -1.0;
114 /* otherwise call gsl */
115 return gsl_sf_erf(*x);
118double FC_FUNC_(oct_legendre_sphplm, OCT_LEGENDRE_SPHPLM)(const fint *l,
119 const int *m,
120 const double *x) {
121 return gsl_sf_legendre_sphPlm(*l, *m, *x);
124/* generalized Laguerre polynomials */
125double FC_FUNC_(oct_sf_laguerre_n, OCT_SF_LAGUERRE_N)(const int *n,
126 const double *a,
127 const double *x) {
128 return gsl_sf_laguerre_n(*n, *a, *x);
131/* Combinations */
133void FC_FUNC_(oct_combination_init, OCT_COMBINATION_INIT)(gsl_combination **c,
134 const fint *n,
135 const fint *k) {
136 *c = gsl_combination_calloc(*n, *k);
137}
138
139void FC_FUNC_(oct_combination_next, OCT_COMBINATION_NEXT)(gsl_combination **c,
140 fint *next) {
141 *next = gsl_combination_next(((gsl_combination *)(*c)));
143
144void FC_FUNC_(oct_get_combination, OCT_GET_COMBINATION)(gsl_combination **c,
145 fint *comb) {
146 int i;
147 for (i = 0; i < ((gsl_combination *)(*c))->k; i++) {
148 comb[i] = (fint)((gsl_combination *)(*c))->data[i];
149 }
152void FC_FUNC_(oct_combination_end, OCT_COMBINATION_END)(gsl_combination **c) {
153 gsl_combination_free(((gsl_combination *)(*c)));
156/* Random Number Generation */
157void FC_FUNC_(oct_ran_init, OCT_RAN_INIT)(gsl_rng **r) {
158 gsl_rng_env_setup();
159 *r = gsl_rng_alloc(gsl_rng_default);
160}
161
162void FC_FUNC_(oct_ran_end, OCT_RAN_END)(gsl_rng **r) { gsl_rng_free(*r); }
163
164/* Random Number Distributions */
165double FC_FUNC_(oct_ran_gaussian, OCT_RAN_GAUSSIAN)(const gsl_rng **r,
166 const double *sigma) {
167 return gsl_ran_gaussian(*r, *sigma);
168}
170double FC_FUNC_(oct_ran_flat, OCT_RAN_FLAT)(const gsl_rng **r, const double *a,
171 const double *b) {
172 return gsl_ran_flat(*r, *a, *b);
175void FC_FUNC_(oct_strerror, OCT_STRERROR)(const fint *err,
176 STR_F_TYPE res STR_ARG1) {
177 const char *c;
178
179 c = gsl_strerror(*err);
180 TO_F_STR1(c, res);
181}
int fint
Definition: fortran_types.h:14
real(real64), dimension(:,:), allocatable sigma
S_E matrix.
Definition: pcm.F90:267
ptrdiff_t l
Definition: operate_inc.c:12
ptrdiff_t i
Definition: operate_inc.c:12