TBCI Numerical high perf. C++ Library  2.8.0
specfun.cpp
Go to the documentation of this file.
1 
5 // $Id: specfun.cpp,v 1.3.2.11 2019/05/28 11:13:02 garloff Exp $
6 
7 //#define NO_NS
8 
9 #include "tbci/specfun.h"
10 #include <stdio.h>
11 #include "tbci/specfun/prototypes2.h"
12 
14 
16  const cplx<double>& z)
17 {
18  printf ("????\n");
19  doublecomplex sin_erg;
20  doublecomplex pow_erg;
21  doublecomplex temp;
23 
24  temp.r = pi*b.real(); temp.i = pi*b.imag();
25  z_sin(&sin_erg, &temp);
26 
27  temp.r = -b.real() + 1.0; temp.i = -b.imag();
28  pow_zz(&pow_erg, (doublecomplex*)(&z), &temp);
29 
30  res = pi/cplx<double>(sin_erg.r, sin_erg.i)
31  *(HypergeometricM(a,b,z)/(gamma(a-b+1)*gamma(b))
32  -cplx<double>(pow_erg.r, pow_erg.i)
33  *HypergeometricM(a-b+1,-b+2,z)/(gamma(a)*gamma(-b+2))
34  );
35  return res;
36 }
37 
38 
40 {
41  floatcomplex arg = {(float)z.real(), (float)z.imag()};
42  floatcomplex erg;
43  cgamma_(&erg, &arg);
44 
45  return cplx<double>(erg.r, erg.i);
46 }
47 
48 
50  const cplx<double>& z)
51 {
52  doublecomplex erg = {0, 0};
53  static const LONG_ int lnchf = 0; //Log-Ausgabe (1)
54  static const LONG_ int ip = 10; //Anzahl der genauen Stellen
55 
56  conhyp_(&erg, (doublecomplex*)&a, (doublecomplex*)&b,
57  (doublecomplex*)&z, &lnchf, &ip);
58 
59  return cplx<double>(erg.r, erg.i);
60 }
61 
62 cplx<double> besselh1(double order, const cplx<double>& z)
63 {
64  double zr, zi, fnu;
65  LONG_ int kode, m, n;
66  double cyr[1], cyi[1];
67  LONG_ int nz, ierr;
68 
69  zr = z.real(); zi = z.imag();
70 
71  kode = 1; // Unskaliert
72  m = 1; // Art H1
73  n = 1; // Laenge der Sequenz
74  fnu = MATH__ fabs(order);
75  zbesh_(&zr, &zi, &fnu, &kode, &m, &n, cyr, cyi, &nz, &ierr);
76  if (ierr || nz)
77  fprintf (stderr, "Error computing Hankel function\n");
78  if (order >= 0.0)
79  return cplx<double> (cyr[0], cyi[0]);
80  //else
81  static const cplx<double> I(0,1);
82  return cplx<double>(cyr[0], cyi[0])* MATH__ exp(pi*fnu*I);
83 }
84 
85 
86 cplx<double> besselh2(double order, const cplx<double>& z)
87 {
88  double zr, zi, fnu;
89  LONG_ int kode, m, n;
90  double cyr[1], cyi[1];
91  LONG_ int nz, ierr;
92 
93  zr = z.real(); zi = z.imag();
94 
95  kode = 1; // Unskaliert
96  m = 2; // Art H2
97  n = 1; // Laenge der Sequenz
98  fnu = MATH__ fabs(order);
99  zbesh_(&zr, &zi, &fnu, &kode, &m, &n, cyr, cyi, &nz, &ierr);
100  if (ierr || nz)
101  fprintf (stderr, "Error computing Hankel function\n");
102  if (order >= 0.0)
103  return cplx<double> (cyr[0], cyi[0]);
104  //else
105  static const cplx<double> I(0,1);
106  return cplx<double> (cyr[0], cyi[0]) * MATH__ exp(-pi*fnu*I);
107 }
108 
109 
110 cplx<double> besselj(double order, const cplx<double>& z)
111 {
112  double zr, zi, fnu;
113  LONG_ int kode, n;
114  double cyr[1], cyi[1];
115  LONG_ int nz, ierr;
116 
117  zr = z.real(); zi = z.imag();
118 
119  kode = 1; // Unskaliert
120  n = 1; // Laenge der Sequenz
121  fnu = MATH__ fabs(order);
122  zbesj_(&zr, &zi, &fnu, &kode, &n, cyr, cyi, &nz, &ierr);
123  if (ierr || nz)
124  fprintf (stderr, "Error computing Hankel function\n");
125  if (order >= 0.0)
126  return cplx<double> (cyr[0], cyi[0]);
127  //else
128  //J(-FNU,Z) = J(FNU,Z)*COS(PI*FNU) - Y(FNU,Z)*SIN(PI*FNU)
129  return cplx<double>(cyr[0], cyi[0]) * MATH__ cos(pi*fnu)
130  - bessely(fnu, z)* MATH__ sin(pi*fnu);
131 }
132 
133 
134 
135 cplx<double> besseli(double order, const cplx<double>& z)
136 {
137  double zr, zi, fnu;
138  LONG_ int kode, n;
139  double cyr[1], cyi[1];
140  LONG_ int nz, ierr;
141 
142  zr = z.real(); zi = z.imag();
143 
144  kode = 1; // Unskaliert
145  n = 1; // Laenge der Sequenz
146  fnu = MATH__ fabs(order);
147  zbesi_(&zr, &zi, &fnu, &kode, &n, cyr, cyi, &nz, &ierr);
148  if (ierr || nz)
149  fprintf (stderr, "Error computing Hankel function\n");
150  if (order >= 0.0)
151  return cplx<double> (cyr[0], cyi[0]);
152  //else
153  //I(-FNU,Z) = I(FNU,Z) + (2/PI)*SIN(PI*FNU)*K(FNU,Z)
154  return cplx<double> (cyr[0], cyi[0])
155  + 2.0/pi*MATH__ sin(pi*fnu)*besselk(fnu,z);
156 }
157 
158 
159 
160 cplx<double> besselk(double order, const cplx<double>& z)
161 {
162  double zr, zi, fnu;
163  LONG_ int kode, n;
164  double cyr[1], cyi[1];
165  LONG_ int nz, ierr;
166 
167  zr = z.real(); zi = z.imag();
168 
169  kode = 1; // Unskaliert
170  n = 1; // Laenge der Sequenz
171  fnu = MATH__ fabs(order);
172  zbesk_(&zr, &zi, &fnu, &kode, &n, cyr, cyi, &nz, &ierr);
173  if (ierr || nz)
174  fprintf (stderr, "Error computing Hankel function\n");
175  if (order >= 0.0)
176  return cplx<double> (cyr[0], cyi[0]);
177  //else
178  //o.k.
179  return cplx<double> (cyr[0], cyi[0]);
180 }
181 
182 
183 cplx<double> bessely(double order, const cplx<double>& z)
184 {
185  double zr, zi, fnu;
186  LONG_ int kode, n;
187  double cyr[1], cyi[1];
188  double wcyr[1], wcyi[1];
189  LONG_ int nz, ierr;
190 
191  zr = z.real(); zi = z.imag();
192 
193  kode = 1; // Unskaliert
194  n = 1; // Laenge der Sequenz
195  fnu = MATH__ fabs(order);
196  zbesy_(&zr, &zi, &fnu, &kode, &n, cyr, cyi, &nz, wcyr, wcyi, &ierr);
197  if (ierr || nz)
198  fprintf (stderr, "Error computing Hankel function\n");
199  if (order >= 0.0)
200  return cplx<double> (cyr[0], cyi[0]);
201  //else
202  //Y(-FNU,Z) = Y(FNU,Z)*COS(PI*FNU) + J(FNU,Z)*SIN(PI*FNU)
203  return cplx<double> (cyr[0], cyi[0])*MATH__ cos(pi*fnu)
204  + besselj(fnu,z)*MATH__ sin(pi*fnu);
205 }
206 
T imag() const
Definition: cplx.h:107
real i
Definition: f2c.h:33
void cgamma_(complex *ret_val, complex *z)
Definition: TOMS404.C:30
int zbesi_(double *zr, double *zi, double *fnu, int *kode, int *n, double *cyr, double *cyi, int *nz, int *ierr)
Definition: zbesi.c:24
Our own complex class.
Definition: cplx.h:48
TBCI__ cplx< T > cos(const TBCI__ cplx< T > &z)
Definition: cplx.h:786
double fabs(const int a)
Definition: basics.h:1215
int zbesk_(double *zr, double *zi, double *fnu, int *kode, int *n, double *cyr, double *cyi, int *nz, int *ierr)
Definition: zbesk.c:25
#define NAMESPACE_TBCI
Definition: basics.h:317
cplx< double > besselk(double order, const cplx< double > &z)
Definition: specfun.cpp:160
void pow_zz(doublecomplex *r, const doublecomplex *a, const doublecomplex *b)
int zbesh_(double *zr, double *zi, double *fnu, int *kode, int *m, int *n, double *cyr, double *cyi, int *nz, int *ierr)
Definition: zbesh.c:132
cplx< double > gamma(const cplx< double > &z)
Definition: specfun.cpp:39
T real() const
Definition: cplx.h:106
int zbesy_(double *zr, double *zi, double *fnu, int *kode, int *n, double *cyr, double *cyi, int *nz, double *wcyr, double *wcyi, int *ierr)
Definition: zbesy.c:23
doublereal r
Definition: f2c.h:34
cplx< double > besseli(double order, const cplx< double > &z)
Definition: specfun.cpp:135
const Vector< T > Vector< T > Vector< T > Vector< T > Vector< T > & z
Definition: LM_fit.h:172
T arg(const TBCI__ cplx< T > &c)
Definition: cplx.h:690
int zbesj_(double *zr, double *zi, double *fnu, int *kode, int *n, double *cyr, double *cyi, int *nz, int *ierr)
Definition: zbesj.c:24
const double pi
Definition: constants.h:38
cplx< double > HypergeometricM(const cplx< double > &a, const cplx< double > &b, const cplx< double > &z)
Definition: specfun.cpp:49
NAMESPACE_TBCI cplx< double > besselh1(double order, const cplx< double > &z)
Definition: specfun.cpp:62
TBCI__ cplx< T > sin(const TBCI__ cplx< T > &z)
Definition: cplx.h:776
F_TMatrix< T > b
Definition: f_matrix.h:736
doublereal i
Definition: f2c.h:34
TBCI__ cplx< T > exp(const TBCI__ cplx< T > &z)
Definition: cplx.h:756
cplx< double > bessely(double order, const cplx< double > &z)
Definition: specfun.cpp:183
void conhyp_(doublecomplex *ret_val, const doublecomplex *a, const doublecomplex *b, const doublecomplex *z, const int *lnchf, const int *ip)
Definition: TOMS_707.C:121
#define LONG_
Definition: prototypes2.h:14
#define NAMESPACE_END
Definition: basics.h:323
cplx< double > besselj(double order, const cplx< double > &z)
Definition: specfun.cpp:110
cplx< double > HypergeometricU(const cplx< double > &a, const cplx< double > &b, const cplx< double > &z)
Definition: specfun.cpp:15
const Vector< T > Vector< T > Vector< T > Vector< T > Vector< T > Vector< T > Vector< T > Vector< T > long int res
Definition: LM_fit.h:199
Definition: f2c.h:33
cplx< double > besselh2(double order, const cplx< double > &z)
Definition: specfun.cpp:86
const unsigned TMatrix< T > const Matrix< T > * a
void z_sin(doublecomplex *r, const doublecomplex *z)
#define MATH__
Definition: basics.h:339
real r
Definition: f2c.h:33