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#include <string.h>
27
28#include <gsl/gsl_combination.h>
29#include <gsl/gsl_randist.h>
30#include <gsl/gsl_rng.h>
31#include <gsl/gsl_sf.h>
32
33#include <fortran_types.h>
34
35/* Numerical threshold for oct_bessel_k0 and oct_bessel_k1 */
36#define BESSEL_K_THRES 1.0e2
37
38/* ---------------------- Interface to GSL functions ------------------------ */
39
40/* Special Functions */
41double FC_FUNC_(oct_incomplete_gamma, OCT_INCOMPLETE_GAMMA)(const double *a,
42 const double *x) {
43 return gsl_sf_gamma_inc_Q(*a, *x);
44}
45
46double FC_FUNC_(oct_sph_bessel, OCT_SPH_BESSEL)(const fint *l,
47 const double *x) {
48 if (*l == 0) {
49 return gsl_sf_bessel_j0(*x);
50 } else if (*l == 1) {
51 return gsl_sf_bessel_j1(*x);
52 } else if (*l == 2) {
53 return gsl_sf_bessel_j2(*x);
54 } else {
55 return gsl_sf_bessel_jl(*l, *x);
56 }
57}
58
59double FC_FUNC_(oct_bessel, OCT_BESSEL)(const fint *n, const double *x) {
60 return gsl_sf_bessel_Jn(*n, *x);
61}
62
63double FC_FUNC_(oct_bessel_in, OCT_BESSEL_IN)(const fint *n, const double *x) {
64 return gsl_sf_bessel_In(*n, *x);
67double FC_FUNC_(oct_bessel_j0, OCT_BESSEL_J0)(const double *x) {
68 return gsl_sf_bessel_J0(*x);
69}
70
71double FC_FUNC_(oct_bessel_j1, OCT_BESSEL_J1)(const double *x) {
72 return gsl_sf_bessel_J1(*x);
76double FC_FUNC_(oct_bessel_k0, OCT_BESSEL_K0)(const double *x) {
77 if (*x > BESSEL_K_THRES) {
78 return 0.0e0;
79 } else {
80 return gsl_sf_bessel_K0(*x);
81 }
82}
83
84double FC_FUNC_(oct_bessel_k1, OCT_BESSEL_K1)(const double *x) {
85 if (*x > BESSEL_K_THRES) {
86 return 0.0e0;
87 } else {
88 return gsl_sf_bessel_K1(*x);
89 }
92double FC_FUNC_(oct_legendre_sphplm, OCT_LEGENDRE_SPHPLM)(const fint *l,
93 const int *m,
94 const double *x) {
95 return gsl_sf_legendre_sphPlm(*l, *m, *x);
96}
98/* generalized Laguerre polynomials */
99double FC_FUNC_(oct_sf_laguerre_n, OCT_SF_LAGUERRE_N)(const int *n,
100 const double *a,
101 const double *x) {
102 return gsl_sf_laguerre_n(*n, *a, *x);
103}
104
105/* Combinations */
107void FC_FUNC_(oct_combination_init, OCT_COMBINATION_INIT)(gsl_combination **c,
108 const fint *n,
109 const fint *k) {
110 *c = gsl_combination_calloc(*n, *k);
111}
112
113void FC_FUNC_(oct_combination_next, OCT_COMBINATION_NEXT)(gsl_combination **c,
114 fint *next) {
115 *next = gsl_combination_next(((gsl_combination *)(*c)));
118void FC_FUNC_(oct_get_combination, OCT_GET_COMBINATION)(gsl_combination **c,
119 fint *comb) {
120 int i;
121 for (i = 0; i < ((gsl_combination *)(*c))->k; i++) {
122 comb[i] = (fint)((gsl_combination *)(*c))->data[i];
125
126void FC_FUNC_(oct_combination_end, OCT_COMBINATION_END)(gsl_combination **c) {
127 gsl_combination_free(((gsl_combination *)(*c)));
130/* Random Number Generation */
131void FC_FUNC_(oct_ran_init, OCT_RAN_INIT)(gsl_rng **r) {
132 gsl_rng_env_setup();
133 *r = gsl_rng_alloc(gsl_rng_default);
134}
136void FC_FUNC_(oct_ran_end, OCT_RAN_END)(gsl_rng **r) { gsl_rng_free(*r); }
137
138/* Random Number Distributions */
139double FC_FUNC_(oct_ran_gaussian, OCT_RAN_GAUSSIAN)(const gsl_rng **r,
140 const double *sigma) {
141 return gsl_ran_gaussian(*r, *sigma);
143
144double FC_FUNC_(oct_ran_flat, OCT_RAN_FLAT)(const gsl_rng **r, const double *a,
145 const double *b) {
146 return gsl_ran_flat(*r, *a, *b);
147}
148
149void oct_strerror(const fint *err, char * res) {
150 const char *c;
151 c = gsl_strerror(*err);
152 strcpy(res, c);
153}
int fint
Definition: fortran_types.h:14
real(real64), dimension(:,:), allocatable sigma
S_E matrix.
Definition: pcm.F90:270
double FC_FUNC_(oct_incomplete_gamma, OCT_INCOMPLETE_GAMMA) const
Definition: oct_gsl_f.c:12269
int fint
Definition: oct_gsl_f.c:12256
void oct_strerror(const fint *err, char *res)
Definition: oct_gsl_f.c:12377
ptrdiff_t l
ptrdiff_t i