Blame view

doc/fvn.tex 70.5 KB
9158e74d6   daniau   git-svn-id: https...
1
2
3
4
5
6
7
  %\documentclass[a4paper,10pt]{article}
  \documentclass[a4paper,english]{article}
  
  \usepackage[utf8]{inputenc}
  \usepackage{a4wide}
  \usepackage{eurosym}
  \usepackage{url}
b40f93a8d   daniau   git-svn-id: https...
8
  \usepackage[colorlinks=true,hyperfigures=true]{hyperref}
9158e74d6   daniau   git-svn-id: https...
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
  %\usepackage{aeguill}
  
  \usepackage{graphicx}
  \usepackage{babel}
  \makeatother
  
  
  %opening
  \title{FVN Documentation}
  \author{William Daniau}
  
  
  \begin{document}
  
  \maketitle
  
  %\begin{abstract}
  
  %\end{abstract}
  \tableofcontents
b93026039   daniau   git-svn-id: https...
29
  \section{Whatis fvn,license}
9158e74d6   daniau   git-svn-id: https...
30
  \subsection{Whatis fvn}
2919a9e2d   daniau   git-svn-id: https...
31
  fvn is a Fortran95 mathematical library with several modules. It provides various usefull subroutine covering linear algebra, numerical integration, least square polynomial, spline interpolation, zero finding, special functions etc.
9158e74d6   daniau   git-svn-id: https...
32

b93026039   daniau   git-svn-id: https...
33
  Most of the work for linear algebra is done by interfacing Lapack \url{http://www.netlib.org/lapack} which means that Lapack and Blas \url{http://www.netlib.org/blas} must be available on your system for linking fvn. If you use an AMD microprocessor, the good idea is to use ACML ( AMD Core Math Library \url{http://developer.amd.com/acml.jsp} which contains an optimized Blas/Lapack.
9158e74d6   daniau   git-svn-id: https...
34

b93026039   daniau   git-svn-id: https...
35
  fvn include some integrated libraries : integration tasks uses a slightly modified version of Quadpack \url{http://www.netlib.org/quadpack}, the fnlib library \url{http://www.netlib.org/fn} is used for special functions and sparse system resolution uses SuiteSparse \url{http://www.cise.ufl.edu/research/sparse/SuiteSparse/}.
9158e74d6   daniau   git-svn-id: https...
36

b93026039   daniau   git-svn-id: https...
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
  This library has been initially written for the use of the ``Acoustic and microsonic'' group leaded by Sylvain Ballandras in the Time and Frequency Department of institute Femto-ST \url{http://www.femto-st.fr/}.
  
  \subsection{License}
      Your use or distribution of fvn or any modified version of
      fvn implies that you agree to this License.
  
      This library is free software; you can redistribute it and/or
      modify it under the terms of the GNU Lesser General Public
      License as published by the Free Software Foundation; either
      version 3 of the License, or (at your option) any later version.
  
      This library 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
      Lesser General Public License for more details.
  
      You should have received a copy of the GNU Lesser General Public
      License along with this library; if not, write to the Free Software
      Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301
      USA
  
      Permission is hereby granted to use or copy this program under the
      terms of the GNU LGPL, provided that the Copyright, this License,
      and the Availability of the original version is retained on all copies.
      User documentation of any code that uses this code or any modified
      version of this code must cite the Copyright, this License, the
      Availability note, and "Used by permission." Permission to modify
      the code and to distribute modified code is granted, provided the
      Copyright, this License, and the Availability note are retained,
      and a notice that the code was modified is included.
9158e74d6   daniau   git-svn-id: https...
67

26c292bb6   wdaniau   Contribution de S...
68
69
70
71
  \subsubsection*{Authors and Contributors}
  \begin{itemize}
   \item William Daniau (william.daniau@femto-st.fr)
   \item Sylvain Ballandras (sylvain.ballandras@femto-st.fr)
8ba5c9c78   wdaniau   1) Updated docume...
72
   \item Christian Waterkeyn (christian.waterkeyn@epcos.com)
26c292bb6   wdaniau   Contribution de S...
73
  \end{itemize}
9158e74d6   daniau   git-svn-id: https...
74

9158e74d6   daniau   git-svn-id: https...
75
76
77
  \section{Naming scheme and convention}
  The naming scheme of the routines is as follow :
  \begin{verbatim}
972644943   daniau   git-svn-id: https...
78
      fvn_*_name()
9158e74d6   daniau   git-svn-id: https...
79
  \end{verbatim} 
972644943   daniau   git-svn-id: https...
80
  where * can be s,d,c or z. 
9158e74d6   daniau   git-svn-id: https...
81
82
83
84
85
86
87
  \begin{itemize}
   \item s is for single precision real (real,real*4,real(4),real(kind=4))
   \item d for double precision real (double precision,real*8,real(8),real(kind=8))
   \item c for single precision complex (complex,complex*8,complex(4),complex(kind=4))
   \item z for double precision complex (double complex,complex*16,complex(8),complex(kind=8))
  \end{itemize}
  In the following description of subroutines parameters, input parameters are followed by (in), output parameters by (out) and parameters which are used as input and modified by the subroutine are followed by (inout).
38581db0c   daniau   git-svn-id: https...
88
  For each routine, there is a generic interface (simply remove \verb'_*' in the name), so using the specific routine is not mandatory.
2919a9e2d   daniau   git-svn-id: https...
89
90
  
  There's a general module called ``\verb'fvn''' that include all fvn submodules. So whatever part of the library is used in the program a ``\verb'use fvn''' will be sufficient instead of specifying the specific module.
8ba5c9c78   wdaniau   1) Updated docume...
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
  \section{kind specification}
  In \verb'fvn' the following definitions are made (this is done in module \verb'fvn_common') to ease portability:
  
  \begin{verbatim}
  integer, parameter :: ip_kind = kind(1)
  integer, parameter :: sp_kind = kind(1.0E0)
  integer, parameter :: dp_kind = kind(1.0D0)
  \end{verbatim}
  
  Although not mandatory, it is a good idea to use these definitions when programming with fvn, that is for example : 
  \begin{verbatim}
  real(kind=sp_kind) :: x 
  real(kind=dp_kind) :: y
  complex(kind=dp_kind) :: z
  \end{verbatim}
  instead of
  \begin{verbatim}
  real :: x
  double precision :: y
  complex(8) :: z
  \end{verbatim}
9158e74d6   daniau   git-svn-id: https...
112
113
114
115
  \section{Linear algebra}
  The linear algebra routines of fvn are an interface to lapack, which make it easier to use.
  \subsection{Matrix inversion}
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
116
  Module : use fvn_linear
38581db0c   daniau   git-svn-id: https...
117
  call fvn_matinv(d,a,inva,status)
9158e74d6   daniau   git-svn-id: https...
118
119
120
  \end{verbatim}
  \begin{itemize}
   \item d (in) is an integer equal to the matrix rank
38581db0c   daniau   git-svn-id: https...
121
122
123
   \item a (in) is a real or complex matrix. It will remain untouched.
   \item inva (out) is a real or complex matrix which contain the inverse of a at the end of the routine
   \item status (out) is an optional integer equal to zero if something went wrong
9158e74d6   daniau   git-svn-id: https...
124
125
126
127
128
  \end{itemize}
  
  \subsubsection*{Example}
  \begin{verbatim}
  program inv
2919a9e2d   daniau   git-svn-id: https...
129
   use fvn_linear
9158e74d6   daniau   git-svn-id: https...
130
131
132
   implicit none
  
   real(8),dimension(3,3) :: a,inva
9158e74d6   daniau   git-svn-id: https...
133
134
135
  
   call random_number(a)
   a=a*100
38581db0c   daniau   git-svn-id: https...
136
   call fvn_matinv(3,a,inva)
9158e74d6   daniau   git-svn-id: https...
137
138
139
140
141
142
143
144
145
146
147
148
   write (*,*) a
   write (*,*)
   write (*,*) inva
   write (*,*)
   write (*,*) matmul(a,inva)
  end program
  \end{verbatim}
  
  
  
  \subsection{Matrix determinants}
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
149
  Module : use fvn_linear
38581db0c   daniau   git-svn-id: https...
150
  det=fvn_det(d,a,status)
9158e74d6   daniau   git-svn-id: https...
151
152
153
  \end{verbatim}
  \begin{itemize}
   \item d (in) is an integer equal to the matrix rank
38581db0c   daniau   git-svn-id: https...
154
155
   \item a (in) is a real or complex matrix. It will remain untouched.
   \item status (out) is an optional integer equal to zero if something went wrong
9158e74d6   daniau   git-svn-id: https...
156
157
158
159
160
  \end{itemize}
  
  \subsubsection*{Example}
  \begin{verbatim}
  program det
2919a9e2d   daniau   git-svn-id: https...
161
   use fvn_linear
9158e74d6   daniau   git-svn-id: https...
162
163
164
165
166
167
168
169
   implicit none
  
   real(8),dimension(3,3) :: a
   real(8) :: deta
   integer :: status
  
   call random_number(a)
   a=a*100
38581db0c   daniau   git-svn-id: https...
170
   deta=fvn_det(3,a,status)
9158e74d6   daniau   git-svn-id: https...
171
172
173
174
175
176
177
178
179
180
181
   write (*,*) a
   write (*,*)
   write (*,*) "Det = ",deta
  end program
  
  \end{verbatim}
  
  
  
  \subsection{Matrix condition}
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
182
  Module : use fvn_linear
38581db0c   daniau   git-svn-id: https...
183
  call fvn_matcon(d,a,rcond,status)
9158e74d6   daniau   git-svn-id: https...
184
185
186
  \end{verbatim}
  \begin{itemize}
   \item d (in) is an integer equal to the matrix rank
38581db0c   daniau   git-svn-id: https...
187
   \item a (in) is a real or complex matrix. It will remain untouched.
9158e74d6   daniau   git-svn-id: https...
188
   \item rcond (out) is a real of same kind as matrix a, it will contain the reciprocal condition number of the matrix
38581db0c   daniau   git-svn-id: https...
189
   \item status (out) is an optional integer equal to zero if something went wrong
9158e74d6   daniau   git-svn-id: https...
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
  \end{itemize}
  
  The reciprocal condition number is evaluated using the 1-norm and is define as in equation \ref{rconddef}
  \begin{equation}
   R = \frac{1}{norm(A)*norm(invA)}
   \label{rconddef}
  \end{equation}
  
  The 1-norm itself is defined as the maximum value of the columns absolute values (modulus for complex) sum as in equation \ref{l1norm}
  \begin{equation}
   L1 = max_j ( \sum_i{\mid A(i,j)\mid}  )
   \label{l1norm}
  \end{equation}
  
  \subsubsection*{Example}
  \begin{verbatim}
  program cond
2919a9e2d   daniau   git-svn-id: https...
207
   use fvn_linear
9158e74d6   daniau   git-svn-id: https...
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
   implicit none
  
   real(8),dimension(3,3) :: a
   real(8) :: rcond
   integer :: status
  
   call random_number(a)
   a=a*100
  
   call fvn_d_matcon(3,a,rcond,status)
   write (*,*) a
   write (*,*)
   write (*,*) "Cond = ",rcond
  end program
  
  \end{verbatim}
  
  
  \subsection{Eigenvalues/Eigenvectors}
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
228
  Module : use fvn_linear
af246bca4   wdaniau   +) Sorting of eig...
229
  call fvn_matev(d,a,evala,eveca,status,sortval)
9158e74d6   daniau   git-svn-id: https...
230
231
232
  \end{verbatim}
  \begin{itemize}
   \item d (in) is an integer equal to the matrix rank
38581db0c   daniau   git-svn-id: https...
233
   \item a (in) is a real or complex matrix. It will remain untouched.
9158e74d6   daniau   git-svn-id: https...
234
235
   \item evala (out) is a complex array of same kind as a. It contains the eigenvalues of matrix a
   \item eveca (out) is a complex matrix of same kind as a. Its columns are the eigenvectors of matrix a : eveca(:,j)=jth eigenvector associated with eigenvalue evala(j).
38581db0c   daniau   git-svn-id: https...
236
   \item status (out) is an optional integer equal to zero if something went wrong
af246bca4   wdaniau   +) Sorting of eig...
237
   \item sortval (in) is an optional logical, if it is true the eigenvalues (and eigenvectors) are sorted in decreasing order of eigenvalues's modulus.
9158e74d6   daniau   git-svn-id: https...
238
  \end{itemize}
2919a9e2d   daniau   git-svn-id: https...
239

9158e74d6   daniau   git-svn-id: https...
240
241
242
  \subsubsection*{Example}
  \begin{verbatim}
  program eigen
2919a9e2d   daniau   git-svn-id: https...
243
   use fvn_linear
9158e74d6   daniau   git-svn-id: https...
244
245
246
247
248
249
250
251
252
   implicit none
  
   real(8),dimension(3,3) :: a
   complex(8),dimension(3) :: evala
   complex(8),dimension(3,3) :: eveca
   integer :: status,i,j
  
   call random_number(a)
   a=a*100
38581db0c   daniau   git-svn-id: https...
253
   call fvn_matev(3,a,evala,eveca,status)
9158e74d6   daniau   git-svn-id: https...
254
255
256
257
258
259
260
261
262
263
264
265
266
267
   write (*,*) a
   write (*,*)
   do i=1,3
      write(*,*) "Eigenvalue ",i,evala(i)
      write(*,*) "Associated Eigenvector :"
      do j=1,3
          write(*,*) eveca(j,i)
      end do
      write(*,*)
   end do
  
  end program
  
  \end{verbatim}
ba393a269   wdaniau   1) updated docume...
268
269
  \subsection{Sparse matrix}
  By interfacing Tim Davis's SuiteSparse from university of Florida \url{http://www.cise.ufl.edu/research/sparse/SuiteSparse/} which is a reference for this kind of problems, fvn provides simple subroutines for solving linear sparse systems and calculating determinants.
9158e74d6   daniau   git-svn-id: https...
270

ba393a269   wdaniau   1) updated docume...
271
  \subsubsection{Sparse solving}
9158e74d6   daniau   git-svn-id: https...
272
273
274
  The provided routines solves the equation $Ax=B$ where A is sparse and given in its triplet form.
  
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
275
  Module : fvn_sparse
ba393a269   wdaniau   1) updated docume...
276
  call fvn_sparse_solve(n,nz,T,Ti,Tj,B,x,status,det)
9158e74d6   daniau   git-svn-id: https...
277
278
  \end{verbatim}
  \begin{itemize}
ba393a269   wdaniau   1) updated docume...
279
280
281
282
283
284
285
286
287
   \item For this family of subroutines the two letters (zl,zi,dl,di) of the specific interface name decribe the arguments's type. z is for complex(8), d for real(8), l for integer(8) and i for integer(4)
   \item \texttt{n} (in) is an integer equal to the matrix rank
   \item \texttt{nz} (in) is an integer equal to the number of non-zero elements
   \item \texttt{T(nz)} (in) is a complex/real array containing the non-zero elements
   \item \textt{Ti(nz)},\texttt{Tj(nz)} (in) are the indexes of the corresponding element of \texttt{T} in the original matrix.
   \item \texttt{B(n)} (in) is a complex/real array containing the second member of the equation.
   \item \texttt{x(n)} (out) is a complex/real array containing the solution
   \item \texttt{status} (out) is an integer which contain non-zero is something went wrong
   \item \texttt{det} (out), is an optional real(8) array of dimension 2 for dl and di specific interface (real systems) and dimension 3 for zl and zi interface (complex systems)
9158e74d6   daniau   git-svn-id: https...
288
  \end{itemize}
2919a9e2d   daniau   git-svn-id: https...
289

9158e74d6   daniau   git-svn-id: https...
290
291
292
  \subsubsection*{Example}
  \begin{verbatim}
  program test_sparse
2919a9e2d   daniau   git-svn-id: https...
293
   use fvn_sparse
9158e74d6   daniau   git-svn-id: https...
294
295
296
297
298
299
300
301
302
303
304
305
306
307
   implicit none
  
   integer(8), parameter :: nz=12
   integer(8), parameter :: n=5
   complex(8),dimension(nz) :: A
   integer(8),dimension(nz) :: Ti,Tj
   complex(8),dimension(n) :: B,x
   integer(8) :: status
  
   A = (/ (2.,0.),(3.,0.),(3.,0.),(-1.,0.),(4.,0.),(4.,0.),(-3.,0.),&
            (1.,0.),(2.,0.),(2.,0.),(6.,0.),(1.,0.) /)
   B = (/ (8.,0.), (45.,0.), (-3.,0.), (3.,0.), (19.,0.)/)
   Ti = (/ 1,2,1,3,5,2,3,4,5,3,2,5 /)
   Tj = (/ 1,1,2,2,2,3,3,3,3,4,5,5 /)
38581db0c   daniau   git-svn-id: https...
308
309
310
   !specific routine that will be used here
   !call fvn_zl_sparse_solve(n,nz,A,Ti,Tj,B,x,status)
   call fvn_sparse_solve(n,nz,A,Ti,Tj,B,x,status)
9158e74d6   daniau   git-svn-id: https...
311
312
313
314
315
316
   write(*,*) x
  
  end program
  
  
  program test_sparse
2919a9e2d   daniau   git-svn-id: https...
317
  use fvn_sparse
9158e74d6   daniau   git-svn-id: https...
318
319
320
321
322
323
324
325
326
327
328
329
330
  implicit none
  
  integer(4), parameter :: nz=12
  integer(4), parameter :: n=5
  real(8),dimension(nz) :: A
  integer(4),dimension(nz) :: Ti,Tj
  real(8),dimension(n) :: B,x
  integer(4) :: status
  
  A = (/ 2.,3.,3.,-1.,4.,4.,-3.,1.,2.,2.,6.,1. /)
  B = (/ 8., 45., -3., 3., 19./)
  Ti = (/ 1,2,1,3,5,2,3,4,5,3,2,5 /)
  Tj = (/ 1,1,2,2,2,3,3,3,3,4,5,5 /)
38581db0c   daniau   git-svn-id: https...
331
332
333
  !specific routine that will be used here
  !call fvn_di_sparse_solve(n,nz,A,Ti,Tj,B,x,status)
  call fvn_sparse_solve(n,nz,A,Ti,Tj,B,x,status)
9158e74d6   daniau   git-svn-id: https...
334
335
336
  write(*,*) x
  
  end program
ba393a269   wdaniau   1) updated docume...
337
338
339
  \end{verbatim}
  
  If optional parameter \texttt{det} is given, the routine will also calculates the matrix determinant and returns it on a mantissa + exponent form, that is the actual determinant will be $det(1).10^{det(2)}$ for real problems and $(det(1)+i.det(2)).10^{det(3)}$ for complex problems. This is given in this form as the determinant can be considerably higher/lower than the biggest/lowest usable double precision real. There's an example of how to use this in following paragraph.
9158e74d6   daniau   git-svn-id: https...
340

ba393a269   wdaniau   1) updated docume...
341
  \subsubsection{Sparse determinant}
9158e74d6   daniau   git-svn-id: https...
342

ba393a269   wdaniau   1) updated docume...
343
344
345
346
  The provided subroutines calculates the determinant of a matrix given in its triplet form.
  \begin{verbatim}
  Module : fvn_sparse
  call fvn_sparse_det(n,nz,T,Ti,Tj,det,status)
9158e74d6   daniau   git-svn-id: https...
347
  \end{verbatim}
ba393a269   wdaniau   1) updated docume...
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
  \begin{itemize}
   \item For this family of subroutines the two letters (zl,zi,dl,di) of the specific interface name decribe the arguments's type. z is for complex(8), d for real(8), l for integer(8) and i for integer(4)
   \item \texttt{n} (in) is an integer equal to the matrix rank
   \item \texttt{nz} (in) is an integer equal to the number of non-zero elements
   \item \texttt{T(nz)} (in) is a complex/real array containing the non-zero elements
   \item \textt{Ti(nz)},\texttt{Tj(nz)} (in) are the indexes of the corresponding element of \texttt{T} in the original matrix.
   \item \texttt{det} (out), a real(8) array of dimension 2 for dl and di specific interface (real systems) and dimension 3 for zl and zi interface (complex systems)
  \item \texttt{status} (out) is an integer which contain non-zero is something went wrong
  \end{itemize}
  
  The matrix determinant is returned on a mantissa + exponent form, that is the actual determinant will be $det(1).10^{det(2)}$ for real problems and $(det(1)+i.det(2)).10^{det(3)}$ for complex problems. This is given in this form as the determinant can be considerably higher/lower than the biggest/lowest usable double precision real.
  
  Here are the possibly returned errors in \texttt{status} parameter :
  \begin{itemize}
   \item \texttt{0} : no errors
   \item \texttt{-1}: out of memory
   \item \texttt{1} : singular matrix
   \item \texttt{2} : determinant underflow, the ``natural'' form of the determinant $det(1).10^{det(2)}$ or $(det(1)+i.det(2)).10^{det(3)}$ will underflow.
   \item \texttt{3} : determinant overflow, the ``natural'' form of the determinant (as above) will overflow
  \end{itemize}
  
  
  And here's an example using this
  
  \begin{verbatim}
  program test_sparse
  use fvn
  implicit none
  integer(kind=sp_kind), parameter :: nz=12
  integer(kind=sp_kind), parameter :: n=5
  complex(kind=dp_kind),dimension(nz) :: A
  complex(kind=dp_kind),dimension(n,n) :: As
  integer(kind=sp_kind),dimension(nz) :: Ti,Tj
  complex(kind=dp_kind),dimension(n) :: B,x
  integer(kind=sp_kind) :: status,i
  real(kind=dp_kind),dimension(3) :: det
  character(len=80) :: fmcmplx
  
  fmcmplx='(5("(",f8.5,",",f8.5,")  "))'
  
  ! Description of the matrix in triplet form
  A = (/ (2.,-1.),(3.,2.),(3.,1.),(-1.,5.),(4.,-7.),(4.,0.),(-3.,-4.),(1.,3.),(2.,0.),(2.,-2.),(6.,4.),(1.,0.) /)
  B = (/ (8.,3.), (45.,1.), (-3.,-2.), (3.,0.), (19.,2.) /)
  Ti = (/ 1,2,1,3,5,2,3,4,5,3,2,5 /)
  Tj = (/ 1,1,2,2,2,3,3,3,3,4,5,5 /)
  
  ! Reconstruction of the matrix in standard form
  As=0.
  do i=1,nz
      As(Ti(i),Tj(i))=A(i)
  end do
  
  write(*,*) "Matrix in standard representation :"
  do i=1,5
      write(*,fmcmplx) As(i,:)
  end do
  write(*,*)
  write(*,*) "Standard determinant : ",fvn_det(5,As)
  write(*,*)
  write(*,*) "Right hand side :"
  write(*,fmcmplx) B
  
  ! can use either specific interface, fvn_zi_sparse_det
  ! either generic one fvn_sparse_det
  call fvn_zi_sparse_det(n,nz,A,Ti,Tj,det,status)
  write(*,*)
  write(*,*) "Sparse Det = ",cmplx(det(1),det(2),kind=dp_kind)*10**det(3)
  ! can use either specific interface fvn_zi_sparse_solve
  ! either generic one fvn_sparse_solve
  ! parameter det is optional
  call fvn_zi_sparse_solve(n,nz,A,Ti,Tj,B,x,status,det)
  write(*,*)
  write(*,*) "Sparse Det as solve option= ",cmplx(det(1),det(2),kind=dp_kind)*10**det(3)
  write(*,*)
  write(*,*) "Solution :"
  write(*,fmcmplx) x
  write(*,*)
  write(*,*) "Product matrix Solution :"
  write(*,fmcmplx) matmul(As,x)
  end program
  
  \end{verbatim}
5b79a897e   daniau   git-svn-id: https...
430
431
  \subsection{Identity matrix}
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
432
   Module : use fvn_linear
81b5d24e1   daniau   git-svn-id: https...
433
   I=fvn_*_ident(n)     (*=s,d,c,z)
5b79a897e   daniau   git-svn-id: https...
434
435
436
437
  \end{verbatim}
  \begin{itemize}
   \item n (in) is an integer equal to the matrix rank
  \end{itemize}
2919a9e2d   daniau   git-svn-id: https...
438

80a3b2e0a   daniau   git-svn-id: https...
439
  This function return the identity matrix of rank n, in the specified type. No generic interface for this one.
5b79a897e   daniau   git-svn-id: https...
440

b93026039   daniau   git-svn-id: https...
441
442
  \subsection{Operators}
  fvn defines some linear operators similar to those defined in IMSL\textregistered (\url{http://www.vni.com/products/imsl/}), that can be used for matrix operations.
2919a9e2d   daniau   git-svn-id: https...
443
444
445
  \begin{verbatim}
   Module : use fvn_linear
  \end{verbatim}
b93026039   daniau   git-svn-id: https...
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
  
  \subsubsection{Unary operators}
  \paragraph{.i.}
  This operator gives the inverse matrix of the argument which must be a square matrix. The status of the operation can be found in the module variable fvn\_status (fourth parameter fvn\_matinv).
  \[
   \textrm{b=.i.a}~\Longleftrightarrow~b=a^{-1}
  \]
  
  \paragraph{.t.}
  This operator gives the transpose matrix of the argument.
  \[
   \textrm{b=.t.a}~\Longleftrightarrow~b= {}^t \! a
  \]
  
  
  \paragraph{.h.}
  This operator gives the conjugate transpose matrix of the argument (also called Hermitian transpose or adjoint matrix).
  \[
   \textrm{b=.h.a}~\Longleftrightarrow~b=a^*=\overline{({}^t\!a)}={}^t \overline{a}
  \]
  
  
  \subsubsection{Binary operators}
  \paragraph{.x.}
  This operator gives the matrix product of the two operands.
  \[
   \textrm{c=a.x.b}~\Longleftrightarrow~c=ab
  \]
  
  
  \paragraph{.ix.}
  This operator gives the matrix product of the inverse of the first operand and the second one. The status of the inversion can be found in the module variable fvn\_status (fourth parameter fvn\_matinv).
  \[
   \textrm{c=a.ix.b}~\Longleftrightarrow~c=a^{-1}b
  \]
  
  \paragraph{.xi.}
  This operator gives the matrix product of the first operand and the inverse of the second one. The status of the inversion can be found in the module variable fvn\_status (fourth parameter fvn\_matinv).
  \[
   \textrm{c=a.xi.b}~\Longleftrightarrow~c=ab^{-1}
  \]
  
  \paragraph{.tx.}
  This operator gives the matrix product of the transpose matrix of the first operand and the second one.
  \[
   \textrm{c=a.tx.b}~\Longleftrightarrow~c={}^t\!ab
  \]
  
  
  \paragraph{.xt.}
  This operator gives the matrix product of the first operand and the transpose matrix of the second one.
  \[
   \textrm{c=a.xt.b}~\Longleftrightarrow~c=a{}^t\!b
  \]
  
  \paragraph{.hx.}
  This operator gives the matrix product of the conjugate transpose matrix of the first operand and the second one.
  \[
   \textrm{c=a.hx.b}~\Longleftrightarrow~c=a^*b
  \]
  
  \paragraph{.xh.}
  This operator gives the matrix product of thee first operand and the conjugate transpose matrix of the second one.
  \[
   \textrm{c=a.xh.b}~\Longleftrightarrow~c=ab^*
  \]
9158e74d6   daniau   git-svn-id: https...
512
513
514
515
516
517
518
519
  
  
  \section{Interpolation}
  
  \subsection{Quadratic Interpolation}
  fvn provide function for interpolating values of a tabulated function of 1, 2 or 3 variables, for both single and double precision.
  \subsubsection{One variable function}
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
520
   Module : use fvn_interpol
38581db0c   daniau   git-svn-id: https...
521
   value=fvn_quad_interpol(x,n,xdata,ydata)
9158e74d6   daniau   git-svn-id: https...
522
523
524
525
526
527
528
529
530
531
532
533
  \end{verbatim}
  \begin{itemize}
   \item x is the real where we want to evaluate the function
   \item n is the number of tabulated values
   \item xdata(n) contains the tabulated coordinates
   \item ydata(n) contains the tabulated function values ydata(i)=y(xdata(i))
  \end{itemize}
  xdata must be strictly increasingly ordered.
  x must be within the range of xdata to actually perform an interpolation, otherwise the resulting value is an extrapolation
  \paragraph*{Example}
  \begin{verbatim}
  program inter1d
2919a9e2d   daniau   git-svn-id: https...
534
  use fvn_interpol
9158e74d6   daniau   git-svn-id: https...
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
  implicit none
  
  integer(kind=4),parameter :: ndata=33
  integer(kind=4) :: i,nout
  real(kind=8) :: f,fdata(ndata),h,pi,q,sin,x,xdata(ndata)
  real(kind=8) ::tv
  
  intrinsic sin
  
  f(x)=sin(x)
  
  xdata(1)=0.
  fdata(1)=f(xdata(1))
  h=1./32.
  do i=2,ndata
        xdata(i)=xdata(i-1)+h
        fdata(i)=f(xdata(i))
  end do
  call random_seed()
  call random_number(x)
  
  q=fvn_d_quad_interpol(x,ndata,xdata,fdata)
  
  tv=f(x)
  write(*,*) "x ",x
  write(*,*) "Calculated (real) value :",tv
  write(*,*) "fvn interpolation :",q
  write(*,*) "Relative fvn error :",abs((q-tv)/tv)
  
  end program
  
  \end{verbatim}
  
  
  \subsubsection{Two variables function}
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
571
  Module : use fvn_interpol
38581db0c   daniau   git-svn-id: https...
572
  value=fvn_quad_2d_interpol(x,y,nx,xdata,ny,ydata,zdata)
9158e74d6   daniau   git-svn-id: https...
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
  \end{verbatim}
  \begin{itemize}
   \item x,y are the real coordinates where we want to evaluate the function
   \item nx is the number of tabulated values along x axis
   \item xdata(nx) contains the tabulated x
   \item ny is the number of tabulated values along y axis
   \item ydata(ny) contains the tabulated y
   \item zdata(nx,ny) contains the tabulated function values zdata(i,j)=z(xdata(i),ydata(j))
  \end{itemize}
  xdata and ydata must be strictly increasingly ordered.
  (x,y) must be within the range of xdata and ydata to actually perform an interpolation, otherwise the resulting value is an extrapolation
  
  \paragraph*{Example}
  
  \begin{verbatim}
  program inter2d
2919a9e2d   daniau   git-svn-id: https...
589
  use fvn_interpol
9158e74d6   daniau   git-svn-id: https...
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
  implicit none
  
  integer(kind=4),parameter  :: nx=21,ny=42
  integer(kind=4) :: i,j
  real(kind=8) :: f,fdata(nx,ny),dble,pi,q,sin,x,xdata(nx),y,ydata(ny)
  real(kind=8) :: tv
  
  intrinsic dble,sin
  
  f(x,y)=sin(x+2.*y)
  do i=1,nx
        xdata(i)=dble(i-1)/dble(nx-1)
  end do
  do i=1,ny
        ydata(i)=dble(i-1)/dble(ny-1)
  end do
  do i=1,nx
        do j=1,ny
              fdata(i,j)=f(xdata(i),ydata(j))
        end do
  end do
  call random_seed()
  call random_number(x)
  call random_number(y)
  
  q=fvn_d_quad_2d_interpol(x,y,nx,xdata,ny,ydata,fdata)
  tv=f(x,y)
  
  write(*,*) "x y",x,y
  write(*,*) "Calculated (real) value :",tv
  write(*,*) "fvn interpolation :",q
  write(*,*) "Relative fvn error :",abs((q-tv)/tv)
  
  end program
  
  \end{verbatim}
  
  
  
  \subsubsection{Three variables function}
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
631
   Module : use fvn_interpol
38581db0c   daniau   git-svn-id: https...
632
  value=fvn_quad_3d_interpol(x,y,z,nx,xdata,ny,ydata,nz,zdata,tdata)
9158e74d6   daniau   git-svn-id: https...
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
  \end{verbatim}
  \begin{itemize}
   \item x,y,z are the real coordinates where we want to evaluate the function
   \item nx is the number of tabulated values along x axis
   \item xdata(nx) contains the tabulated x
   \item ny is the number of tabulated values along y axis
   \item ydata(ny) contains the tabulated y
   \item nz is the number of tabulated values along z axis
   \item zdata(ny) contains the tabulated z
   \item tdata(nx,ny,nz) contains the tabulated function values tdata(i,j,k)=t(xdata(i),ydata(j),zdata(k))
  \end{itemize}
  xdata, ydata and zdata must be strictly increasingly ordered.
  (x,y,z) must be within the range of xdata and ydata to actually perform an interpolation, otherwise the resulting value is an extrapolation
  
  \paragraph*{Example}
  \begin{verbatim}
  program inter3d
2919a9e2d   daniau   git-svn-id: https...
650
  use fvn_interpol
9158e74d6   daniau   git-svn-id: https...
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
  
  implicit none
  
  integer(kind=4),parameter  :: nx=21,ny=42,nz=18
  integer(kind=4) :: i,j,k
  real(kind=8) :: f,fdata(nx,ny,nz),dble,pi,q,sin,x,xdata(nx),y,ydata(ny),z,zdata(nz)
  real(kind=8) :: tv
  
  intrinsic dble,sin
  
  f(x,y,z)=sin(x+2.*y+3.*z)
  do i=1,nx
        xdata(i)=2.*(dble(i-1)/dble(nx-1))
  end do
  do i=1,ny
        ydata(i)=2.*(dble(i-1)/dble(ny-1))
  end do
  do i=1,nz
        zdata(i)=2.*(dble(i-1)/dble(nz-1))
  end do
  do i=1,nx
        do j=1,ny
              do k=1,nz
                    fdata(i,j,k)=f(xdata(i),ydata(j),zdata(k))
              end do
        end do
  end do
  call random_seed()
  call random_number(x)
  call random_number(y)
  call random_number(z)
  
  q=fvn_d_quad_3d_interpol(x,y,z,nx,xdata,ny,ydata,nz,zdata,fdata)
  tv=f(x,y,z)
  
  write(*,*) "x y z",x,y,z
  write(*,*) "Calculated (real) value :",tv
  write(*,*) "fvn interpolation :",q
  write(*,*) "Relative fvn error :",abs((q-tv)/tv)
  
  end program
  
  \end{verbatim}
  
  \subsubsection{Utility procedure}
  fvn provides a simple utility procedure to locate the interval in which a value is located in an increasingly ordered array.
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
698
  Module : use fvn_interpol
38581db0c   daniau   git-svn-id: https...
699
  call fvn_find_interval(x,i,xdata,n)
9158e74d6   daniau   git-svn-id: https...
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
  \end{verbatim}
  \begin{itemize}
   \item x (in) the real value to locate
   \item i (out) the resulting indice
   \item xdata(n) (in) increasingly ordered array
   \item n (in) size of the array
  \end{itemize}
  The resulting integer i is as : $xdata(i) <= x < xdata(i+1)$. If $x < xdata(1)$ then $i=0$ is returned. If $x > xdata(n)$ then $i=n$ is returned. Finally if $x=xdata(n)$ then $i=n-1$ is returned.
  
  
  
  \subsection{Akima spline}
  fvn provides Akima spline interpolation and evaluation for both single and double precision real.
  \subsubsection{Interpolation}
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
715
  Module : use fvn_interpol
38581db0c   daniau   git-svn-id: https...
716
  call fvn_akima(n,x,y,br,co)
9158e74d6   daniau   git-svn-id: https...
717
718
719
720
721
722
723
724
725
726
  \end{verbatim}
  \begin{itemize}
   \item n (in) is an integer equal to the number of points
   \item x(n) (in) ,y(n) (in) are the known couples of coordinates
   \item br (out) on output contains a copy of x
   \item co(4,n) (out) is a real matrix containing the 4 coefficients of the Akima interpolation spline for a given interval.
  \end{itemize}
  
  \subsubsection{Evaluation}
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
727
  Module : use fvn_interpol
38581db0c   daniau   git-svn-id: https...
728
  y=fvn_spline_eval(x,n,br,co)
9158e74d6   daniau   git-svn-id: https...
729
730
731
732
733
734
735
736
737
738
739
740
741
742
  \end{verbatim}
  \begin{itemize}
   \item x (in) is the point where we want to evaluate
   \item n (in) is the number of known points and br(n) (in), co(4,n) (in) \\
  are the outputs of fvn\_x\_akima(n,x,y,br,co) 
  \end{itemize}
  
  \subsubsection{Example}
  In the following example we will use Akima splines to interpolate a sinus function with 30 points between -10 and 10. We then use the evaluation function to calculate the coordinates of 1000 points between -11 and 11, and write a 3 columns file containing : x, calculated sin(x), interpolation evaluation of sin(x).
  
  One can see that the interpolation is very efficient even with only 30 points. Of course as soon as we leave the -10 to 10 interval, the values are extrapolated and thus can lead to very inacurrate values.
  
  \begin{verbatim}
  program akima
2919a9e2d   daniau   git-svn-id: https...
743
   use fvn_interpol
9158e74d6   daniau   git-svn-id: https...
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
   implicit none
  
   integer :: nbpoints,nppoints,i
   real(8),dimension(:),allocatable :: x_d,y_d,breakpoints_d
   real(8),dimension(:,:),allocatable :: coeff_fvn_d
   real(8) :: xstep_d,xp_d,ty_d,fvn_y_d
  
   open(2,file='fvn_akima_double.dat')
   open(3,file='fvn_akima_breakpoints_double.dat')
   nbpoints=30
   allocate(x_d(nbpoints))
   allocate(y_d(nbpoints))
   allocate(breakpoints_d(nbpoints))
   allocate(coeff_fvn_d(4,nbpoints))
  
   xstep_d=20./dfloat(nbpoints)
   do i=1,nbpoints
      x_d(i)=-10.+dfloat(i)*xstep_d
      y_d(i)=dsin(x_d(i))
      write(3,44) (x_d(i),y_d(i))
   end do
   close(3)
  
   call fvn_d_akima(nbpoints,x_d,y_d,breakpoints_d,coeff_fvn_d)
  
   nppoints=1000 
   xstep_d=22./dfloat(nppoints)
   do i=1,nppoints
      xp_d=-11.+dfloat(i)*xstep_d
      ty_d=dsin(xp_d)
      fvn_y_d=fvn_d_spline_eval(xp_d,nbpoints-1,breakpoints_d,coeff_fvn_d)
      write(2,44) (xp_d,ty_d,fvn_y_d)
   end do
  
   close(2)
  
  44      FORMAT(4(1X,1PE22.14))
  
  end program
  
  \end{verbatim}
  Results are plotted on figure \ref{akima}
  
  \begin{figure}
   \begin{center}
   \includegraphics[width=0.9\textwidth]{akima.pdf}
   % akima.pdf: 504x720 pixel, 72dpi, 17.78x25.40 cm, bb=0 0 504 720
   \caption{Akima Spline Interpolation}
   \label{akima}
  \end{center}
  
  \end{figure}
  
  
  
  \section{Least square polynomial}
80a3b2e0a   daniau   git-svn-id: https...
800
  fvn provide a function to find a least square polynomial of a given degree, for real in single or double precision. It is performed using Lapack subroutine sgels (dgels), which solve this problem.
9158e74d6   daniau   git-svn-id: https...
801
802
  
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
803
  Module : use fvn_linear
38581db0c   daniau   git-svn-id: https...
804
  call fvn_lspoly(np,x,y,deg,coeff,status)
9158e74d6   daniau   git-svn-id: https...
805
806
807
808
809
810
811
812
813
814
815
816
817
  \end{verbatim}
  \begin{itemize}
   \item np (in) is an integer equal to the number of points
   \item x(np) (in),y(np) (in) are the known coordinates
   \item deg (in) is an integer equal to the degree of the desired polynomial, it must be lower than np.
   \item coeff(deg+1) (out) on output contains the polynomial coefficients
   \item status (out) is an integer containing 0 if a problem occured.
  \end{itemize}
  
  \subsection*{Example}
  Here's a simple example : we've got 13 measurement points and we want to find the least square degree 3 polynomial for these points :
  \begin{verbatim}
   program lsp
2919a9e2d   daniau   git-svn-id: https...
818
   use fvn_linear
9158e74d6   daniau   git-svn-id: https...
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
   implicit none
  
   integer,parameter :: npoints=13,deg=3
   integer :: status,i
   real(kind=8) :: xm(npoints),ym(npoints),xstep,xc,yc
   real(kind=8) :: coeff(deg+1)
  
   xm = (/ -3.8,-2.7,-2.2,-1.9,-1.1,-0.7,0.5,1.7,2.,2.8,3.2,3.8,4. /)
   ym = (/ -3.1,-2.,-0.9,0.8,1.8,0.4,2.1,1.8,3.2,2.8,3.9,5.2,7.5  /)
  
   open(2,file='fvn_lsp_double_mesure.dat')
   open(3,file='fvn_lsp_double_poly.dat')
  
   do i=1,npoints
      write(2,44) xm(i),ym(i)
   end do
   close(2)
  
  
   call fvn_d_lspoly(npoints,xm,ym,deg,coeff,status)
  
   xstep=(xm(npoints)-xm(1))/1000.
   do i=1,1000
      xc=xm(1)+(i-1)*xstep
      yc=poly(xc,coeff)
      write(3,44) xc,yc
   end do
   close(3)
  
  44      FORMAT(4(1X,1PE22.14))
  
  contains
  function poly(x,coeff)
      implicit none
      real(8) :: x
      real(8) :: coeff(deg+1)
      real(8) :: poly
      integer :: i
  
      poly=0.
  
      do i=1,deg+1
          poly=poly+coeff(i)*x**(i-1)
      end do
  
  end function
  end program
  \end{verbatim}
  The results are plotted on figure \ref{lsp} .
  
  \begin{figure}
   \begin{center}
   \includegraphics[width=0.9\textwidth]{lsp.pdf}
   \caption{Least Square Polynomial}
   \label{lsp}
   \end{center}
  \end{figure}
  
  
  
  \section{Zero finding}
  fvn provide a routine for finding zeros of a complex function using Muller algorithm (only for double complex type). It is based on a version provided on the web by Hans D Mittelmann \url{http://plato.asu.edu/ftp/other\_software/muller.f}.
  
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
883
   Module : use fvn_misc
38581db0c   daniau   git-svn-id: https...
884
   call fvn_muller(f,eps,eps1,kn,nguess,n,x,itmax,infer,ier)
9158e74d6   daniau   git-svn-id: https...
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
  \end{verbatim}
  \begin{itemize}
   \item f (in) is the complex function (kind=8) for which we search zeros
   \item eps (in) is a real(8) corresponding to the first stopping criterion : let fp(z)=f(z)/p where p = (z-z(1))*(z-z(2))*,,,*(z-z(k-1)) and z(1),...,z(k-1) are previously found roots.  if ((cdabs(f(z)).le.eps) .and. (cdabs(fp(z)).le.eps)), then z is accepted as a root.
   \item eps1 (in) is a real(8) corresponding to the second stopping criterion : a root is accepted if two successive approximations to a given root agree within eps1. Note that if either or both of the stopping criteria are fulfilled, the root is accepted.
   \item kn (in) is an integer equal to the number of known roots, which must be stored in x(1),...,x(kn), prior to entry in the subroutine.
   \item nguess (in) is the number of initial guesses provided. These guesses must be stored in x(kn+1),...,x(kn+nguess). nguess must be set equal to zero if no guesses are provided.
   \item n (in) is an integer equal to the number of new roots to be found.
   \item x (inout) is a complex(8) vector of length kn+n. x(1),...,x(kn) on input must contain any known roots.  x(kn+1),..., x(kn+n) on input may, on user option, contain initial guesses for the n new roots which are to be computed. If the user does not provide an initial guess, zero is used. On output, x(kn+1),...,x(kn+n) contain the approximate roots found by the subroutine.
   \item itmax (in) is an integer equal to the maximum allowable number of iterations per root.
   \item infer (out) is an integer vector of size kn+n. On output infer(j) contains the number of iterations used in finding the j-th root when convergence was achieved.  If convergence was not obtained in itmax iterations, infer(j) will be greater than itmax
   \item ier (out) is an integer used as an error parameter. ier = 33 indicates failure to converge within itmax iterations for at least one of the (n) new roots.
  \end{itemize}
  This subroutine always returns the last approximation for root j in x(j). if the convergence criterion is satisfied, then infer(j) is less than or equal to itmax. if the convergence criterion is not satisified, then infer(j) is set to either itmax+1 or itmax+k, with k greater than 1. infer(j) = itmax+1 indicates that muller did not obtain convergence in the allowed number of iterations. in this case, the user may wish to set itmax to a larger value. infer(j) = itmax+k means that convergence was obtained (on iteration k) for the deflated function fp(z) = f(z)/((z-z(1)...(z-z(j-1))) but failed for f(z). in this case, better initial guesses might help or, it might be necessary to relax the convergence criterion.
  
  \subsection*{Example}
  Example to find the ten roots of $x^{10}-1$
  \begin{verbatim}
   program muller
2919a9e2d   daniau   git-svn-id: https...
904
   use fvn_misc
9158e74d6   daniau   git-svn-id: https...
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
   implicit none
  
   integer :: i,info
   complex(8),dimension(10) :: roots
   integer,dimension(10) :: infer
   complex(8), external :: f
  
   call fvn_z_muller(f,1.d-12,1.d-10,0,0,10,roots,200,infer,info)
  
   write(*,*) "Error code :",info
   do i=1,10
      write(*,*) roots(i),infer(i)
   enddo
   end program
  
   function f(x)
      complex(8) :: x,f
      f=x**10-1
   end function
  
  \end{verbatim}
9158e74d6   daniau   git-svn-id: https...
926
927
928
929
930
931
932
  
  \section{Numerical integration}
  Using an integrated slightly modified version of quadpack \url{http://www.netlib.org/quadpack}, fvn provide adaptative numerical integration (Gauss Kronrod) of real functions of 1 and 2 variables. fvn also provide a function to calculate Gauss-Legendre abscissas and weight, and a simple non adaptative integration subroutine. All routines exists only in fvn for double precision real.
  
  \subsection{Gauss Legendre Abscissas and Weigth}
  This subroutine was inspired by Numerical Recipes routine gauleg.
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
933
  Module : use fvn_integ
38581db0c   daniau   git-svn-id: https...
934
  call fvn_gauss_legendre(n,qx,qw)
9158e74d6   daniau   git-svn-id: https...
935
936
937
938
939
940
941
942
943
944
  \end{verbatim}
  \begin{itemize}
   \item n (in) is an integer equal to the number of Gauss Legendre points
   \item qx (out) is a real(8) vector of length n containing the abscissas.
   \item qw (out) is a real(8) vector of length n containing the weigths.
  \end{itemize}
  This subroutine computes n Gauss-Legendre abscissas and weigths
  
  \subsection{Gauss Legendre Numerical Integration}
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
945
  Module : use fvn_integ
38581db0c   daniau   git-svn-id: https...
946
  call fvn_gl_integ(f,a,b,n,res)
9158e74d6   daniau   git-svn-id: https...
947
948
949
950
951
952
953
954
955
956
957
958
959
  \end{verbatim}
  \begin{itemize}
   \item f (in) is a real(8) function to integrate
   \item a (in) and b (in) are real(8) respectively lower and higher bound of integration
   \item n (in) is an integer equal to the number of Gauss Legendre points to use
   \item res (out) is a real(8) containing the result
  \end{itemize}
  This function is a simple Gauss Legendre integration subroutine, which evaluate the integral of function f as in equation \ref{intsple} using n Gauss-Legendre pairs.
  
  \subsection{Gauss Kronrod Adaptative Integration}
  This kind of numerical integration is an iterative procedure which try to achieve a given precision.
  \subsubsection{Numerical integration of a one variable function}
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
960
  Module : use fvn_integ
38581db0c   daniau   git-svn-id: https...
961
  call fvn_integ_1_gk(f,a,b,epsabs,epsrel,key,res,abserr,ier,limit)
9158e74d6   daniau   git-svn-id: https...
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
  \end{verbatim}
  This routine evaluate the integral of function f as in equation \ref{intsple}
  \begin{itemize}
   \item f (in) is an external real(8) function of one variable
   \item a (in) and b (in) are real(8) respectively lower an higher bound of integration
   \item epsabs (in) and epsrel (in) are real(8) respectively desired absolute and relative error
   \item key (in) is an integer between 1 and 6 correspondind to the Gauss-Kronrod rule to use :
      \begin{itemize}
          \item 1 : 7 - 15 points
          \item 2 : 10 - 21 points
          \item 3 : 15 - 31 points
          \item 4 : 20 - 41 points
          \item 5 : 25 - 51 points
          \item 6 : 30 - 61 points
      \end{itemize}
   \item res (out) is a real(8) containing the estimation of the integration.
   \item abserr (out) is a real(8) equal to the estimated absolute error
   \item ier (out) is an integer used as an error flag
      \begin{itemize}
          \item 0 : no error
          \item 1 : maximum number of subdivisions allowed has been achieved. one can allow more subdivisions by increasing the value of limit (and taking the according dimension adjustments into account). however, if this yield no improvement it is advised to analyze the integrand in order to determine the integration difficulaties. If the position of a local difficulty can be determined (i.e.singularity, discontinuity within the interval) one will probably gain from splitting up the interval at this point and calling the integrator on the subranges. If possible, an appropriate special-purpose integrator should be used which is designed for handling the type of difficulty involved.
          \item 2 : the occurrence of roundoff error is detected, which prevents the requested tolerance from being achieved.
          \item 3 : extremely bad integrand behaviour occurs at some points of the integration interval.
          \item 6 : the input is invalid, because (epsabs.le.0 and epsrel.lt.max(50*rel.mach.acc.,0.5d-28)) or limit.lt.1 or lenw.lt.limit*4. result, abserr, neval, last are set to zero. Except when lenw is invalid, iwork(1), work(limit*2+1) and work(limit*3+1) are set to zero, work(1) is set to a and work(limit+1) to b.
      \end{itemize}
93f1019c0   daniau   git-svn-id: https...
987
   \item limit (in) is an optional integer equal to maximum number of subintervals in the partition of the given integration interval (a,b). If the parameter is not present a default value of 500 will be used.
9158e74d6   daniau   git-svn-id: https...
988
989
990
991
992
993
994
995
996
997
998
999
  \end{itemize}
  
  \begin{equation}
   \int_a^b f(x)~dx
   \label{intsple}
  \end{equation}
  
  
  
  
  \subsubsection{Numerical integration of a two variable function}
  \begin{verbatim}
2919a9e2d   daniau   git-svn-id: https...
1000
  Module : use fvn_integ
38581db0c   daniau   git-svn-id: https...
1001
  call fvn_integ_2_gk(f,a,b,g,h,epsabs,epsrel,key,res,abserr,ier,limit)
9158e74d6   daniau   git-svn-id: https...
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
  \end{verbatim}
  This function evaluate the integral of a function f(x,y) as defined in equation \ref{intdble}. The parameters of same name as in the previous paragraph have exactly the same function and behaviour thus only what differs is decribed here
  \begin{itemize}
   \item a (in) and b (in) are real(8) corresponding respectively to lower and higher bound of integration for the x variable.
   \item g(x) (in) and h(x) (in) are external functions describing the lower and higher bound of integration for the y variable as a function of x.
  \end{itemize}
  
  \begin{equation}
   \int_a^b \int_{g(x)}^{h(x)} f(x,y)~dy~dx
   \label{intdble}
  \end{equation}
  
  \subsubsection*{Example}
  \begin{verbatim}
  program integ
2919a9e2d   daniau   git-svn-id: https...
1017
   use fvn_integ
9158e74d6   daniau   git-svn-id: https...
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
   implicit none
  
   real(8), external :: f1,f2,g,h
   real(8) :: a,b,epsabs,epsrel,abserr,res
   integer :: key,ier
  
   a=0.
   b=1.
   epsabs=1d-8
   epsrel=1d-8
   key=2
   call fvn_d_integ_1_gk(f1,a,b,epsabs,epsrel,key,res,abserr,ier,500)
   write(*,*) "Integration of x*x between 0 and 1 : "
   write(*,*) res
  
   call fvn_d_integ_2_gk(f2,a,b,g,h,epsabs,epsrel,key,res,abserr,ier,500)
   write(*,*) "Integration of x*y between 0 and 1 on both x and y : "
   write(*,*) res
   
  
  end program
  
  function f1(x)
   implicit none
      real(8) :: x,f1
      f1=x*x
  end function
  
  function f2(x,y)
   implicit none
      real(8) :: x,y,f2
      f2=x*y
  end function
  
  function g(x)
   implicit none
      real(8) :: x,g
      g=0.
  end function
  
  function h(x)
   implicit none
      real(8) :: x,h
      h=1.
  end function
  \end{verbatim}
5b79a897e   daniau   git-svn-id: https...
1064
  \section{Special functions}
80a3b2e0a   daniau   git-svn-id: https...
1065
  Specials functions are available in fvn by using an implementation of fnlib \url{http://www.netlib.org/fn} with some additions. This can be used separatly from the rest of fvn by using the module \verb'fvn_fnlib' and linking the library \verb'libfvn_fnlib.a' . The module provides a generic interfaces to all the routines. Specific names of the routines are given in the description.
38581db0c   daniau   git-svn-id: https...
1066

2919a9e2d   daniau   git-svn-id: https...
1067
1068
1069
  \begin{verbatim}
   Module : use fvn_fnlib
  \end{verbatim}
38581db0c   daniau   git-svn-id: https...
1070
  \paragraph{Important Note}
2919a9e2d   daniau   git-svn-id: https...
1071
  Due to the addition of fnlib to fvn, some functions that were in fvn and are redondant are now removed from fvn, so update your code now and replace them with the fnlib version. These are listed here after :
38581db0c   daniau   git-svn-id: https...
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
  \begin{itemize}
   \item \verb'fvn_z_acos' replaced by \verb'acos'
   \item \verb'fvn_z_asin' replaced by \verb'asin'
   \item \verb'fvn_d_asinh' replaced by \verb'asinh'
   \item \verb'fvn_d_acosh' replaced by \verb'acosh'
   \item \verb'fvn_s_csevl' replaced by \verb'csevl'
   \item \verb'fvn_d_csevl' replaced by \verb'csevl'
   \item \verb'fvn_d_factorial' replaced by \verb'fac'
   \item \verb'fvn_d_lngamma' replaced by \verb'alngam'
  \end{itemize}
  
  
  \subsection{Elementary functions}
  \subsubsection{carg}
  \begin{verbatim}
   carg(z)
  \end{verbatim}
  \begin{itemize}
   \item z (in) is a complex
  \end{itemize}
  This function evaluates the argument of the complex z. That is $\theta$ for $z=\rho e^{i\theta}$.
  
  Specific interfaces : \verb'carg,zarg'
  
  
  \subsubsection{cbrt}
  \begin{verbatim}
   cbrt(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates the cubic root of the argument x.
  
  Specific interfaces : \verb'cbrt,dcbrt,ccbrt,zcbrt'
  
  
  \subsubsection{exprl}
  \begin{verbatim}
   exprl(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates ${e^x-1}\over x$.
  
  Specific interfaces : \verb'exprel,dexprl,cexprl,zexprl'
  
  \subsubsection{log10}
  \begin{verbatim}
   log10(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function is an extension of the intrinsic function log10 to complex arguments.
  
  Specific interfaces : \verb'clog10,zlog10'
  
  
  \subsubsection{alnrel}
  \begin{verbatim}
   alnrel(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates $ln(1+x)$.
  
  Specific interfaces : \verb'alnrel,dlnrel,clnrel,zlnrel'
  
  
  \subsection{Trigonometry}
  \subsubsection{tan}
  \begin{verbatim}
   tan(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates the tangent of the argument. It is an extension of the intrinsic function tan to complex arguments.
  
  Specific interfaces : \verb'ctan,ztan'
  
  
  \subsubsection{cot}
  \begin{verbatim}
   cot(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluate the cotangent of the argument.
  
  Specific interfaces : \verb'cot,dcot,ccot,zcot'
  
  \subsubsection{sindg}
  \begin{verbatim}
   sindg(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluate the sinus of the argument expressed in degrees.
  
  Specific interfaces : \verb'sindg,dsindg'
  
  
  \subsubsection{cosdg}
  \begin{verbatim}
   cosdg(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluate the cosinus of the argument expressed in degrees.
  
  Specific interfaces : \verb'cosdg,dcosdg'
  
  
  \subsubsection{asin}
  \begin{verbatim}
   asin(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates the arc sine of the argument. It is an extension of the intrinsic function asin to complex arguments.
  
  Specific interfaces : \verb'casin,zasin'
  
  \subsubsection{acos}
  \begin{verbatim}
   acos(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates the arc cosine of the argument. It is an extension of the intrinsic function acos to complex arguments.
  
  Specific interfaces : \verb'cacos,zacos'
  
  
  \subsubsection{atan}
  \begin{verbatim}
   atan(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates the arc tangent of the argument. It is an extension of the intrinsic function atan to complex arguments.
  
  Specific interfaces : \verb'catan,zatan'
  
  \subsubsection{atan2}
  \begin{verbatim}
   atan2(x,y)
  \end{verbatim}
  \begin{itemize}
   \item x,y are real or complex
  \end{itemize}
  This function evaluates the arc tangent of $x \over y$. It is an extension of the intrinsic function atan2 to complex arguments.
  
  Specific interfaces : \verb'catan2,zatan2'
  
  \subsubsection{sinh}
  \begin{verbatim}
   sinh(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates the hyperbolic sine of the argument. It is an extension of the intrinsic function sinh to complex arguments.
  
  Specific interfaces : \verb'csinh,zsinh'
  
  
  \subsubsection{cosh}
  \begin{verbatim}
   cosh(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates the hyperbolic cosine of the argument. It is an extension of the intrinsic function cosh to complex arguments.
  
  Specific interfaces : \verb'ccosh,zcosh'
  
  \subsubsection{tanh}
  \begin{verbatim}
   tanh(x)
  \end{verbatim}
  This function evaluates the hyperbolic tangent of the argument. It is an extension of the intrinsic function tanh to complex arguments.
  
  Specific interfaces : \verb'ctanh,ztanh'
  
  \subsubsection{asinh}
  \begin{verbatim}
   asinh(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates the arc hyperbolic sine of the argument.
  
  Specific interfaces : \verb'asinh,dasinh,casinh,zasinh'
  
  \subsubsection{acosh}
  \begin{verbatim}
   acosh(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates the arc hyperbolic cosine of the argument.
  
  Specific interfaces : \verb'acosh,dacosh,cacosh,zacosh'
  
  \subsubsection{atanh}
  \begin{verbatim}
   atanh(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates the arc hyperbolic tangent of the argument.
  
  Specific interfaces : \verb'atanh,datanh,catanh,zatanh'
  
  \subsection{Exponential Integral and related}
  \subsubsection{ei}
  \begin{verbatim}
   ei(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates the exponential integral for argument greater then 0 and the Cauchy principal value for argument less than 0. It is define by equation \ref{ei} for $x 
  eq 0$.
  \begin{equation}
  \label{ei}
   ei(x)= - \int _{-x} ^\infty {e^{-t}\over t}dt
  \end{equation}
  
  Specific interfaces : \verb'ei,dei'
  
  
  \subsubsection{e1}
  \begin{verbatim}
   e1(x)
  \end{verbatim}
  \begin{itemize}
26c292bb6   wdaniau   Contribution de S...
1324
   \item x is a real or complex
38581db0c   daniau   git-svn-id: https...
1325
  \end{itemize}
26c292bb6   wdaniau   Contribution de S...
1326
1327
  For a real argument, this function evaluates the exponential integral for argument greater than 0 and the Cauchy principal value for argument less than 0. It is define by equation \ref{e1} for $x 
  eq 0$.
38581db0c   daniau   git-svn-id: https...
1328
1329
1330
1331
  \begin{equation}
  \label{e1}
   e1(x)= \int _{x} ^\infty {e^{-t}\over t}dt
  \end{equation}
26c292bb6   wdaniau   Contribution de S...
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
  For a complex argument, the notation in equation \ref{e1cplx} is used (Abramowitz and Stegun, p.228 \url{http://www.math.ucla.edu/~cbm/aands/page_228.htm}):
  \begin{equation}
  \label{e1cplx}
   e1(z) = \int _{z} ^\infty {e^{-t}\over t}dt \textrm{~with~} \left|  arg(z) \right| < \pi
  \end{equation}
  For positive values of real part of $z$, this can be written as in equation \ref{e1cplx_pos} :
  \begin{equation}
  \label{e1cplx_pos}
   e1(z)= \int _{1} ^\infty {e^{-tz}\over t}dt \textrm{~with~} Re(z) > 0
  \end{equation}
  
   
  
  Specific interfaces : \verb'e1,de1,ze1'
38581db0c   daniau   git-svn-id: https...
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
  
  \subsubsection{ali}
  \begin{verbatim}
   ali(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates the logarithm integral. it is define by equation \ref{ali} for $x > 0$ and $x 
  eq 1$.
  \begin{equation}
   \label{ali}
   ali(x)= - \int _0 ^x {dt \over ln(x)}
  \end{equation}
  
  Specific interfaces : \verb'ali,dli'
  
  \subsubsection{si}
  \begin{verbatim}
   si(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates the sine integral defined by equation \ref{si}.
  \begin{equation}
   \label{si}
   si(x)= \int _0 ^x {sin(t) \over t }dt
  \end{equation}
  
  Specific interfaces : \verb'si,dsi'
  
  
  \subsubsection{ci}
  \begin{verbatim}
   ci(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates the cosine integral defined by equation \ref{ci} where $\gamma \approx 0.57721566$ represent Euler's constant.
  \begin{equation}
   \label{ci}
   ci(x)= \gamma + ln(x) + \int _0 ^x {{1-cos(t)} \over t} dt
  \end{equation}
  
  Specific interfaces : \verb'ci,dci'
  
  \subsubsection{cin}
5b79a897e   daniau   git-svn-id: https...
1395
  \begin{verbatim}
38581db0c   daniau   git-svn-id: https...
1396
   cin(x)
5b79a897e   daniau   git-svn-id: https...
1397
1398
  \end{verbatim}
  \begin{itemize}
38581db0c   daniau   git-svn-id: https...
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
   \item x is a real
  \end{itemize}
  This function evaluates the cosine integral alternate definition given by equation \ref{cin}.
  \begin{equation}
   \label{cin}
   cin(x)= \int _0 ^x {{1-cos(t)} \over t} dt
  \end{equation}
  
  Specific interface : \verb'cin,dcin'
  
  \subsubsection{shi}
  \begin{equation}
   shi(x)
  \end{equation}
  \begin{itemize}
   \item x is a real
5b79a897e   daniau   git-svn-id: https...
1415
  \end{itemize}
38581db0c   daniau   git-svn-id: https...
1416
1417
1418
1419
1420
1421
1422
  This function evaluates the hyperbolic sine integral defined by equation \ref{shi}.
  \begin{equation}
   \label{shi}
   shi(x) = \int _0 ^x {sinh(t) \over t}dt
  \end{equation}
  
  Specific interfaces : \verb'shi,dshi'
9158e74d6   daniau   git-svn-id: https...
1423

38581db0c   daniau   git-svn-id: https...
1424
  \subsubsection{chi}
5b79a897e   daniau   git-svn-id: https...
1425
  \begin{verbatim}
38581db0c   daniau   git-svn-id: https...
1426
   chi(x)
5b79a897e   daniau   git-svn-id: https...
1427
1428
  \end{verbatim}
  \begin{itemize}
38581db0c   daniau   git-svn-id: https...
1429
   \item x is a real
5b79a897e   daniau   git-svn-id: https...
1430
  \end{itemize}
38581db0c   daniau   git-svn-id: https...
1431
1432
1433
1434
1435
1436
1437
  This function evaluates the hyperbolic cosine integral defined by equation \ref{chi} where $\gamma \approx 0.57721566$ represent Euler's constant.
  \begin{equation}
   \label{chi}
   chi(x)= \gamma + ln(x) + \int _0 ^x {{cosh(t) -1} \over t}dt
  \end{equation}
  
  Specific interfaces : chi,dchi
5b79a897e   daniau   git-svn-id: https...
1438

38581db0c   daniau   git-svn-id: https...
1439
1440
  
  \subsubsection{cinh}
5b79a897e   daniau   git-svn-id: https...
1441
  \begin{verbatim}
38581db0c   daniau   git-svn-id: https...
1442
   cinh(x)
5b79a897e   daniau   git-svn-id: https...
1443
1444
1445
  \end{verbatim}
  \begin{itemize}
   \item x is a real
38581db0c   daniau   git-svn-id: https...
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
  \end{itemize}
  This function evaluates the hyperbolic cosine integral alternate definition given by equation \ref{cinh}.
  \begin{equation}
   \label{cinh}
   cinh(x) = \int _0 ^x {{cosh(t) -1} \over t}dt
  \end{equation}
  
  Specific interfaces : cinh,dcinh
  
  
  \subsection{Gamma function and related}
  \subsubsection{fac}
  \begin{verbatim}
   fac(n)
81b5d24e1   daniau   git-svn-id: https...
1460
   dfac(n)
38581db0c   daniau   git-svn-id: https...
1461
1462
  \end{verbatim}
  \begin{itemize}
5b79a897e   daniau   git-svn-id: https...
1463
1464
   \item n is an integer
  \end{itemize}
81b5d24e1   daniau   git-svn-id: https...
1465
  This function return $n!$ as a real(4) or real(8) for dfac. There's no generic interface for this one.
38581db0c   daniau   git-svn-id: https...
1466
1467
1468
1469
1470
1471
  
  Specific interfaces : \verb'fac,dfac'
  
  \subsubsection{binom}
  \begin{verbatim}
   binom(n,m)
81b5d24e1   daniau   git-svn-id: https...
1472
   dbinom(n,m)
38581db0c   daniau   git-svn-id: https...
1473
1474
1475
1476
  \end{verbatim}
  \begin{itemize}
   \item n,m are integers
  \end{itemize}
81b5d24e1   daniau   git-svn-id: https...
1477
  This function return the binomial coefficient defined by equation \ref{binom} with $n \geq m \geq 0$. binom returns a real(4), dbinom a real(8). There's no generic interface for this one.
38581db0c   daniau   git-svn-id: https...
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
  \begin{equation}
   \label{binom}
   binom(n,m) = C_n^m = {{n!} \over {m!(n-m)!}}
  \end{equation}
  
  Specific interfaces : \verb'binom,dbinom'
  
  
  \subsubsection{gamma}
  \begin{verbatim}
   gamma(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates $ \Gamma (x) $ defined by equation \ref{gamma}.
  \begin{equation}
   \label{gamma}
   \Gamma (x) = \int _0 ^{\infty} t^{x-1}e^{-t}dt
  \end{equation}
  Note that $n!=\Gamma (n+1)$.
  
  Specific interfaces :\verb'gamma,dgamma,cgamma,zgamm'
  
  \subsubsection{gamr}
  \begin{verbatim}
   gamr(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates the reciprocal gamma function $gamr(x)= {1 \over \Gamma(x)}$
  
  
  \subsubsection{alngam}
  \begin{verbatim}
   alngam(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates $ln(|\Gamma(x)|)$
  
  Specific interfaces : \verb'alngam,dlngam,clngam,zlngam'
  
  
  \subsubsection{algams}
  \begin{verbatim}
   call algams(x,algam,sgngam)
  \end{verbatim}
  \begin{itemize}
   \item x (in) is a real
   \item algam (out) is a real
   \item sgngam (out) is a real
  \end{itemize}
  This subroutine evaluates the logarithm of the absolute value of gamma and the sign of gamma. 
  $algam=ln(|\Gamma(x)|)$ and $sgngam=1.0$ or $-1.0$ according to the sign of $\Gamma(x)$.
  
  Specific interfaces : \verb'algams,dlgams'
  
  \subsubsection{gami}
  \begin{verbatim}
   gami(a,x)
  \end{verbatim}
  \begin{itemize}
   \item x is a positive real
   \item a is a strictly positive real
  \end{itemize}
  This function evaluates the incomplete gamma function defined by equation \ref{gami}.
  \begin{equation}
   \label{gami}
   gami(a,x)=\gamma(a,x)=\int _0 ^x t^{a-1} e^{-t}dt
  \end{equation}
  
  Specific interfaces : \verb'gami,dgami'
  
  \subsubsection{gamic}
  \begin{verbatim}
   gamic(a,x)
  \end{verbatim}
  \begin{itemize}
   \item x is a positive real
   \item a is a real
  \end{itemize}
  This function evaluates the complementary incomplete gamma function defined by equation \ref{gamic}.
  \begin{equation}
   \label{gamic}
   gamic(a,x)=\Gamma(a,x)=\int _x ^\infty t^{a-1} e^{-t}dt
  \end{equation}
  
  Specific interfaces : \verb'gamic,dgamic'
  
  \subsubsection{gamit}
  \begin{verbatim}
   gamit(a,x)
  \end{verbatim}
  \begin{itemize}
   \item x is a positive real
   \item a is a real
  \end{itemize}
  This function evaluates the Tricomi's incomplete gamma function defined by equation \ref{gamit}.
  \begin{equation}
   \label{gamit}
   gamit(a,x)=\gamma^* (a,x)= {{x^{-a}\gamma(a,x)}\over \Gamma(a)}
  \end{equation}
  
  Specific interfaces : \verb'gamit,dgamit'
  
  
  \subsubsection{psi}
  \begin{verbatim}
   psi(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real or complex
  \end{itemize}
  This function evaluates the psi function which is the logarithm derivative of the gamma function as defined in equation \ref{psi}.
  \begin{equation}
   \label{psi}
    psi(x)= \psi(x) = {d\over dx} ln(\Gamma(x))
  \end{equation}
  x must not be zero or a negative integer.
  
  Specific interfaces : \verb'psi,dpsi,cpsi,zpsi'
  
  
  \subsubsection{poch}
  \begin{verbatim}
   poch(a,x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
   \item a is a real
  \end{itemize}
  This function evaluates a generalization of Pochhammer's symbol.
  
  Pochhammer's symbol for n a positive integer is given by equation \ref{poch_int}
  \begin{equation}
   \label{poch_int}
  (a)_n = a(a-1)(a-2)...(a-n+1)
  \end{equation}
  
  The generalization of Pochhammer's symbol is given by equation \ref{poch}
  \begin{equation}
   \label{poch}
   poch(a,x)= (a)_x = {\Gamma(a+x) \over \Gamma(a)}
  \end{equation}
  
  Specific interfaces : \verb'poch,dpoch'
  
  
  \subsubsection{poch1}
  \begin{verbatim}
   poch1(a,x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
   \item a is a real
  \end{itemize}
  This function is defined by equation \ref{poch1}. It is usefull for certains situations, especially when x is small.
  
  \begin{equation}
   \label{poch1}
   poch1(a,x)={{(a)_x-1} \over x}
  \end{equation}
  
  Specific interfaces : \verb'poch1,dpoch1'
  
  \subsubsection{beta}
  \begin{verbatim}
   beta(a,b)
  \end{verbatim}
  \begin{itemize}
   \item a,b are real positive or complex
  \end{itemize}
  This function evaluates $\beta$ function defined by equation \ref{beta}.
  \begin{equation}
   \label{beta}
  beta(a,b)=\beta(a,b)={  {\Gamma(a) \Gamma(b)} \over \Gamma(a+b) }
  \end{equation}
  
  Specific interfaces : \verb'beta,dbeta,cbeta,zbeta'
  
  
  \subsubsection{albeta}
  \begin{verbatim}
   albeta(a,b)
  \end{verbatim}
  \begin{itemize}
   \item a,b are real positive or complex
  \end{itemize}
  This function evaluates the natural logarithm of beta function : $ln(\beta(a,b))$
  
  Specific interfaces : \verb'albeta,dlbeta,clbeta,zlbeta'
  
  \subsubsection{betai}
  \begin{verbatim}
   betai(x,pin,qin)
  \end{verbatim}
  \begin{itemize}
   \item x is a real in [0,1]
   \item pin and qin are strictly positive real
  \end{itemize}
  This function evaluates the incomplete beta function ratio, that is the probability that a random variable from a beta distribution having parameters pin and qin will be less than or equal to x. It is defined by equation \ref{betai}.
  
  \begin{equation}
   \label{betai}
  betai(x,pin,qin)=I_x(pin,qin)={1 \over \beta(pin,qin)} \int _0 ^x t^{pin-1}(1-t)^{qin-1}dt
  \end{equation}
  
  Specific interfaces : \verb'betai,dbetai'
  
  \subsection{Error function and related}
  \subsubsection{erf}
  \begin{verbatim}
   erf(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates the error function defined by equation \ref{erf}.
  \begin{equation}
   \label{erf}
   erf(x)={2\over \sqrt{ \pi}} \int _0 ^x e^{-t^2}dt
  \end{equation}
  
  Specific interfaces : \verb'erf,derf'
  
  \subsubsection{erfc}
  \begin{verbatim}
   erfc(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates the complimentary error function defined by equation \ref{erfc}.
  \begin{equation}
   \label{erfc}
   erfc(x)={2\over \sqrt{ \pi}} \int _x ^\infty e^{-t^2}dt
  \end{equation}
  
  Specific interfaces : \verb'erfc,derfc'
  
  
  \subsubsection{daws}
  \begin{verbatim}
   daws(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates Dawson's function defined by equation \ref{daws}.
  \begin{equation}
   \label{daws}
  daws(x)=e^{-x^2} \int _0 ^x e^{t^2}dt
  \end{equation}
  
  Specific interfaces : \verb'daws,ddaws'
  
  \subsection{Bessel functions and related}
  \subsubsection{bsj0}
  \begin{verbatim}
   bsj0(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates Bessel function of the first kind of order 0 defined by equation \ref{bsj0}.
  \begin{equation}
   \label{bsj0}
   bsj0(x)=J_0(x)= {1 \over \pi} \int _0 ^\pi cos(x sin(\theta)) d\theta
  \end{equation}
  
  Specific interfaces : \verb'besj0,dbesj0'
  
  \subsubsection{bsj1}
  \begin{verbatim}
   bsj1(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates Bessel function of the first kind of order 1 defined by equation \ref{bsj1}.
  \begin{equation}
   \label{bsj1}
  bsj1(x)=J_1(x)={1 \over \pi} \int _0 ^\pi cos(x sin(\theta)-\theta) d\theta
  \end{equation}
  
  Specific interfaces : \verb'besj1,dbesj1'
80a3b2e0a   daniau   git-svn-id: https...
1767
1768
1769
1770
1771
1772
1773
1774
1775
1776
1777
1778
1779
1780
1781
1782
1783
1784
1785
1786
  \subsubsection{bsjn}
  \begin{verbatim}
   bsjn(n,x,factor,big)
  \end{verbatim}
  \begin{itemize}
   \item n is an integer
   \item x is a real
   \item factor is an optional integer
   \item big is an optional real
  \end{itemize}
  This function evaluates Bessel function of the first kind of order n (plotted in figure \ref{besseljnfamily}). These functions satisfy the recurrent relation \ref {bsjn}.
  \begin{equation}
   \label{bsjn}
  J_{n+1}(x)={{2n}\over x} J_n(x) - J_{n-1}(x)
  \end{equation}
  This relation is directly used in upward direction to compute $J_n(x)$ for $x>n$. However it is unstable for $x<n$, therefore a Miller's Algorithm is used. The principle of this method is to use the reccurent relation downward from an arbitrary higher than n order with an arbitrary seed and then normalize the solution with \ref{bsjnnorm}
  \begin{equation}
   \label{bsjnnorm}
   1=J_0+2J_2+2J_4+2J_6+...
  \end{equation}
8ba5c9c78   wdaniau   1) Updated docume...
1787
  The optional parameters \verb'factor' and \verb'big' can be used to modify the behaviour of the algorithm. \verb'factor' is used in determining the arbitrary starting order ( an even integer near $n+\sqrt{factor~n}$), the default $factor$ value is 40 for single precision and 150 for double precision. \verb'big' is a real determining the threshold for which anti-overflow counter measures has to be taken, default value is $1.10^{10}$
80a3b2e0a   daniau   git-svn-id: https...
1788
1789
1790
1791
1792
1793
1794
1795
1796
1797
1798
1799
1800
  
  By convenience, the routine accept $n=0$ and $n=1$, in that cases a call to $bsj0(x)$ or $bsj1(x)$ is actually performed.
  
  \begin{figure}
   \begin{center}
   \includegraphics[width=0.9\textwidth]{besselj-mono.pdf}
   \caption{Bessel $J_n$ functions family}
   \label{besseljnfamily}
   \end{center}
  \end{figure}
  
  
  Specific interfaces : \verb'besjn,dbesjn'
8ba5c9c78   wdaniau   1) Updated docume...
1801
1802
1803
1804
1805
1806
1807
1808
1809
1810
1811
1812
1813
1814
1815
  
  \subsubsection{besrj}
  \begin{verbatim}
   call besrj(x,n,b)
  \end{verbatim}
  \begin{itemize}
   \item x (in) is a real
   \item n (in) is an integer
   \item b (out) is a real array of dimension n
  \end{itemize}
  This subroutine evaluates Bessel function of the first kind of order \verb'0' to \verb'n-1' for argument x and return the result in array b, which then contain \verb'b(1)='$J_0(x)$,\verb'b(2)='$J_1(x)$,...,\verb'b(n)='$J_{n-1}(x)$.
  
  The algorithm is different from the one used in \verb'bsjn', the choice between the two depends on accuracy and timing considerations, there are test programs (\verb'test_besrj' and \verb'test_bestime') in the \verb'fvn_test' directory which can help choosing the good one. 
  
  Specific interfaces : \verb'besrj,dbesrj'
38581db0c   daniau   git-svn-id: https...
1816
1817
1818
1819
1820
  \subsubsection{bsy0}
  \begin{verbatim}
   bsy0(x)
  \end{verbatim}
  \begin{itemize}
80a3b2e0a   daniau   git-svn-id: https...
1821
   \item x is a strictly positive real
38581db0c   daniau   git-svn-id: https...
1822
1823
1824
1825
1826
1827
1828
1829
1830
1831
1832
1833
1834
1835
  \end{itemize}
  This function evaluates the Bessel function of the second kind of order 0 defined by equation \ref{bsy0}
  \begin{equation}
   \label{bsy0}
   bsy0(x)=Y_0(x)={1 \over \pi} \int _0 ^\pi sin(x sin(\theta))d\theta -{2 \over \pi} \int _0 ^\infty e^{-x sinh(t)}dt
  \end{equation}
  
  Specific interfaces : \verb'besy0,dbesy0'
  
  \subsubsection{bsy1}
  \begin{verbatim}
   bsy1(x)
  \end{verbatim}
  \begin{itemize}
80a3b2e0a   daniau   git-svn-id: https...
1836
   \item x is a strictly positive real
38581db0c   daniau   git-svn-id: https...
1837
1838
1839
1840
1841
1842
1843
1844
1845
  \end{itemize}
  This function evaluates the Bessel function of the second kind of order 1 defined by equation \ref{bsy1}.
  \begin{equation}
   \label{bsy1}
   bsy1(x)=Y_1(x)=-{1 \over \pi} \int _0 ^\pi sin (\theta - x sin(\theta)) d\theta 
  - {1 \over \pi} \int _0 ^\infty (e^t -e^{-t})e^{-x sinh(t)}dt
  \end{equation}
  
  Specific interfaces : \verb'besy1,dbesy1'
80a3b2e0a   daniau   git-svn-id: https...
1846
1847
1848
1849
1850
1851
1852
1853
1854
1855
1856
1857
1858
1859
1860
1861
1862
1863
1864
1865
1866
1867
1868
1869
1870
1871
1872
1873
  
  \subsubsection{bsyn}
  \begin{verbatim}
   bsyn(n,x)
  \end{verbatim}
  \begin{itemize}
   \item n is an integer
   \item x is a strictly positive real
  \end{itemize}
  This function evaluates the Bessel function of the second kind of order n (plotted in figure \ref{besselynfamily}). These functions satisfy the recurrent relation \ref {bsyn}.
  \begin{equation}
   \label{bsyn}
  Y_{n+1}(x)={{2n}\over x} Y_n(x) - Y_{n-1}(x)
  \end{equation}
  This recurrent relation is directly used in the upward direction to compute $Y_n(x)$.
  
  By convenience, the routine accept $n=0$ and $n=1$, in that cases a call to $bsy0(x)$ or $bsy1(x)$ is actually performed.
  
  \begin{figure}
   \begin{center}
   \includegraphics[width=0.9\textwidth]{bessely-mono.pdf}
   \caption{Bessel $Y_n$ functions family}
   \label{besselynfamily}
   \end{center}
  \end{figure}
  
  
  Specific interfaces : \verb'besyn,dbesyn'
38581db0c   daniau   git-svn-id: https...
1874
1875
1876
1877
1878
1879
1880
1881
1882
1883
1884
1885
1886
1887
1888
1889
1890
1891
1892
1893
1894
1895
1896
1897
1898
1899
1900
1901
1902
  \subsubsection{bsi0}
  \begin{verbatim}
   bsi0(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates the Bessel function of the third kind of order 0 defined by equation \ref{bsi0}.
  \begin{equation}
   \label{bsi0}
  bsi0(x)=I_0(x)={1 \over \pi} \int _0 ^\pi cosh(x cos(\theta))d\theta
  \end{equation}
  
  Specific interfaces : \verb'besi0,dbesi0'
  
  \subsubsection{bsi1}
  \begin{verbatim}
   bsi1(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates the Bessel function of the third kind of order 1 defined by equation \ref{bsi1}.
  \begin{equation}
   \label{bsi1}
  bsi1(x)=I_1(x)={1 \over \pi} \int _0 ^\pi e^{x cos(\theta)} cos(\theta)d\theta
  \end{equation}
  
  Specific interfaces : \verb'besi1,dbesi1'
80a3b2e0a   daniau   git-svn-id: https...
1903
1904
1905
1906
1907
1908
1909
1910
1911
1912
1913
1914
1915
1916
1917
1918
1919
1920
1921
1922
  \subsubsection{bsin}
  \begin{verbatim}
   bsin(n,x,factor,big)
  \end{verbatim}
  \begin{itemize}
   \item n is an integer
   \item x is a real
   \item factor is an optional integer
   \item big is an optional real
  \end{itemize}
  This function evaluates Bessel function of the third kind of order n (plotted in figure \ref{besselinfamily}). These functions satisfy the recuurent relation \ref{bsin}
  \begin{equation}
   \label{bsin}
  I_{n+1}(x)=-{{2n}\over x}I_n(x) + I_{n-1}(x)
  \end{equation}
  This relation is unstable in the upward direction, therefore a Miller's Algorithm is used to evaluate the function. Even if there's a usable normalization relation \ref{bsinnorm}, it is not used in the routine, instead normalization  is done by a simple call to \verb'bsi0(x)'.
  \begin{equation}
   \label{bsinnorm}
   1=I_0-2I_2+2I_4-2I_6+....
  \end{equation}
8ba5c9c78   wdaniau   1) Updated docume...
1923
  The optional parameters \verb'factor' and \verb'big' can be used to modify the behaviour of the algorithm. \verb'factor' is used in determining the arbitrary starting order ( an even integer near $n+\sqrt{factor~n}$), the default $factor$ value is 40 for single precision and 150 for double precision. \verb'big' is a real determining the threshold for which anti-overflow counter measures has to be taken, default value is $1.10^{10}$
80a3b2e0a   daniau   git-svn-id: https...
1924
1925
1926
1927
1928
1929
1930
1931
1932
1933
1934
1935
1936
  
  By convenience, the routine accept $n=0$ and $n=1$, in that cases a call to $bsi0(x)$ or $bsi1(x)$ is actually performed.
  
  \begin{figure}
   \begin{center}
   \includegraphics[width=0.9\textwidth]{besseli-mono.pdf}
   \caption{Bessel $I_n$ functions family}
   \label{besselinfamily}
   \end{center}
  \end{figure}
  
  
  Specific interfaces : \verb'besin,dbesin'
8ba5c9c78   wdaniau   1) Updated docume...
1937
1938
1939
1940
1941
1942
1943
1944
1945
1946
1947
1948
1949
1950
  \subsubsection{besri}
  \begin{verbatim}
   call besri(x,n,b)
  \end{verbatim}
  \begin{itemize}
   \item x (in) is a real
   \item n (in) is an integer
   \item b (out) is a real array of dimension n
  \end{itemize}
  This subroutine evaluates Bessel function of the third kind of order \verb'0' to \verb'n-1' for argument x and return the result in array b, which then contain \verb'b(1)='$I_0(x)$,\verb'b(2)='$I_1(x)$,...,\verb'b(n)='$I_{n-1}(x)$.
  
  The algorithm is different from the one used in \verb'bsin' and tends to be more accurate, the choice between the two depends on accuracy and timing considerations, there are test programs (\verb'test_besri' and \verb'test_bestime') in the \verb'fvn_test' directory which can help choosing the good one. 
  
  Specific interfaces : \verb'besri,dbesri'
38581db0c   daniau   git-svn-id: https...
1951
1952
1953
1954
1955
1956
1957
1958
1959
1960
1961
1962
1963
1964
1965
1966
1967
1968
1969
1970
1971
1972
1973
1974
1975
1976
1977
1978
1979
  \subsubsection{bsk0}
  \begin{verbatim}
   bsk0(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a strictly positive real
  \end{itemize}
  This function evaluates the modified Bessel function of the second kind of order 0 defined by equation \ref{bsk0}
  \begin{equation}
   \label{bsk0}
  bsk0(x)=K_0(x)=\int _0 ^\infty cos(x sinh(t))dt
  \end{equation}
  
  Specific interfaces : \verb'besk0,dbesk0'
  
  \subsubsection{bsk1}
  \begin{verbatim}
   bsk1(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a strictly positive real
  \end{itemize}
  This function evaluates the modified Bessel function of the second kind of order 1 defined by equation \ref{bsk1}
  \begin{equation}
   \label{bsk1}
  bsk1(x)=K_1(x)=\int _0 ^\infty sin(x sinh(t))sinh(t)dt
  \end{equation}
  
  Specific interfaces : \verb'besk1,dbesk1'
80a3b2e0a   daniau   git-svn-id: https...
1980
1981
1982
1983
1984
1985
1986
1987
1988
1989
1990
1991
1992
1993
1994
1995
1996
1997
1998
1999
2000
2001
2002
2003
2004
2005
2006
2007
  \subsubsection{bskn}
  \begin{verbatim}
   bskn(n,x)
  \end{verbatim}
  \begin{itemize}
   \item n is an integer
   \item x is strictly positive real
  \end{itemize}
  This function evaluates the modified Bessel function of the second kind of order n (plotted in figure \ref{besselknfamily}). These functions satisfy the recurrent relation \ref{bskn}
  
  \begin{equation}
   \label{bskn}
  K_{n+1}(x)={{2n}\over x}K_n(x)+K_{n-1}(x)
  \end{equation}
  This recurrent relation is directly used in the upward direction to compute $K_n(x)$.
  
  By convenience, the routine accept $n=0$ and $n=1$, in that cases a call to $bsk0(x)$ or $bsk1(x)$ is actually performed.
  
  \begin{figure}
   \begin{center}
   \includegraphics[width=0.9\textwidth]{besselk-mono.pdf}
   \caption{Bessel $K_n$ functions family}
   \label{besselknfamily}
   \end{center}
  \end{figure}
  
  
  Specific interface : \verb'beskn,dbeskn'
38581db0c   daniau   git-svn-id: https...
2008
2009
2010
2011
2012
2013
2014
2015
2016
2017
2018
2019
2020
2021
2022
2023
2024
2025
2026
2027
2028
2029
2030
2031
2032
2033
2034
2035
2036
2037
2038
2039
2040
2041
2042
2043
2044
2045
2046
2047
2048
2049
2050
2051
2052
2053
2054
2055
2056
2057
2058
2059
2060
2061
2062
2063
2064
2065
2066
2067
2068
2069
2070
2071
2072
2073
2074
2075
2076
2077
2078
2079
2080
2081
2082
2083
2084
2085
2086
2087
2088
2089
2090
2091
2092
2093
2094
2095
2096
2097
2098
2099
2100
2101
2102
2103
2104
2105
2106
2107
2108
2109
2110
2111
2112
2113
2114
2115
2116
2117
2118
2119
2120
2121
2122
2123
2124
2125
2126
2127
2128
2129
2130
2131
2132
2133
2134
2135
2136
2137
2138
2139
2140
2141
2142
2143
2144
2145
2146
2147
2148
2149
2150
2151
2152
2153
2154
2155
2156
2157
2158
2159
2160
  \subsubsection{bsi0e}
  \begin{verbatim}
   bsi0e(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates $e^{-|x|}I_0(x)$
  
  Specific interfaces : \verb'besi0e,dbsi0e'
  
  \subsubsection{bsi1e}
  \begin{verbatim}
   bsi1e(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates $e^{-|x|}I_1(x)$
  
  Specific interfaces : \verb'besi1e,dbsi1e'
  
  \subsubsection{bsk0e}
  \begin{verbatim}
   bsk0e(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a strictly positive real
  \end{itemize}
  This function evaluates $e^x K_0(x)$
  
  Specific interfaces : \verb'besk0e,dbsk0e'
  
  \subsubsection{bsk1e}
  \begin{verbatim}
   bsk1e(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a strictly positive real
  \end{itemize}
  This function evaluates $e^x K_1(x)$
  
  Specific interfaces : \verb'besk1e,dbsk1e'
  
  \subsubsection{bsks}
  \begin{verbatim}
   call bsks(xnu,x,nin,bk)
  \end{verbatim}
  \begin{itemize}
   \item xnu (in) is a real with $|xnu|<1$. It's the fractional order
   \item x (in) is a real. The value for which the sequence of Bessel functions is to be evaluated.
   \item nin (in) is an integer.
   \item bk (out) is a real vector of length abs(nin), containing the values of the function.
  \end{itemize}
  This subroutine evaluates a sequence of modified Bessel function of the second kind of fractional order.
  
  If nin is positive, on completion $bk(1)=K_
  u(x)$,$bk(2)=K_{
  u+1}(x)$,...,$bk(nin)=K_{
  u+nin-1}(x)$. If nin is negative, on completion $bk(1)=K_
  u(x)$,$bk(2)=K_{
  u-1}(x)$,...,$bk(|nin|)=K_{
  u+nin+1}(x)$.
  
  Specific interfaces : \verb'besks,dbesks'
  
  \subsubsection{bskes}
  \begin{verbatim}
   call bskes(xnu,x,nin,bke)
  \end{verbatim}
  \begin{itemize}
   \item xnu (in) is a real with $|xnu|<1$. It's the fractional order
   \item x (in) is a real. The value for which the sequence of exponentialy scaled Bessel functions is to be evaluated.
   \item nin (in) is an integer. Number of elements in the sequence.
   \item bke (out) is a real vector of length abs(nin), containing the values of the function.
  \end{itemize}
  This subroutine evaluates a sequence of exponentially scaled modified Bessel function of the second kind of fractional order.
  
  If nin is positive, on completion $bk(1)=e^x K_
  u(x)$,$bk(2)=e^x K_{
  u+1}(x)$,...,$bk(nin)=e^x K_{
  u+nin-1}(x)$. If nin is negative, on completion $bk(1)=e^x K_
  u(x)$,$bk(2)=e^x K_{
  u-1}(x)$,...,$bk(|nin|)=e^x K_{
  u+nin+1}(x)$.
  
  Specific interfaces : \verb'beskes,dbskes'
  
  \subsection{Airy function and related}
  \subsubsection{ai}
  \begin{verbatim}
   ai(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates the airy function defined by equation \ref{ai}
  \begin{equation}
   \label{ai}
  Ai(x)={1 \over \pi} \int _0 ^\infty cos(xt+ {1 \over 3}t^3)dt
  \end{equation}
  
  Specific interfaces : \verb'ai,dai'
  
  \subsubsection{bi}
  \begin{verbatim}
   bi(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates the Airy function of the second kind defined by equation \ref{bi}
  \begin{equation}
   \label{bi}
  Bi(x)={1 \over \pi} \int _0 ^\infty e^{xt- {1 \over 3}t^3}dt +
  {1 \over \pi} \int _0 ^\infty sin(xt+ {1 \over 3}t^3)dt
  \end{equation}
  
  Specific interfaces : \verb'bi,dbi'
  
  
  \subsubsection{aid}
  \begin{verbatim}
   aid(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates the derivative of the Airy function, $aid(x)={d \over {dx}}Ai(x)$.
  
  Specific interface : \verb'aid,daid'
  
  \subsubsection{bid}
  \begin{verbatim}
   bid(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates the derivative of the Airy function of the second kind, $bid(x)={d \over {dx}}Bi(x)$.
  
  Specific interfaces : \verb'bid,dbid'
  
  
  \subsubsection{aie}
  \begin{verbatim}
   aie(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  
  This function evaluates the exponentially scaled Airy function defined in equation \ref{aie}.
93f1019c0   daniau   git-svn-id: https...
2161
  %\begin{equation}
80a3b2e0a   daniau   git-svn-id: https...
2162
2163
  % \label{aie}
  %aie(x)=Ai(x) \textrm{~if~}x\leq0 \qquad \qquad aie(x)=e^{{2\over3}x^{3\over2}}Ai(x)\rm{~if~}x>0
93f1019c0   daniau   git-svn-id: https...
2164
  %\end{equation}
80a3b2e0a   daniau   git-svn-id: https...
2165
2166
2167
2168
2169
2170
2171
  \begin{equation}
  \label{aie}
    aie(x)=\left\{\begin{array}{ll}
    Ai(x) & \textrm{~if~}x\leq0 \\
    e^{{2\over3}x^{3\over2}}Ai(x) & \textrm{~if~}x>0
  \end{array}\right.
  \end{equation}
38581db0c   daniau   git-svn-id: https...
2172
2173
2174
2175
2176
2177
2178
2179
2180
2181
2182
2183
  
  Specific interfaces : \verb'aie,daie'
  
  \subsubsection{bie}
  
  \begin{verbatim}
   bie(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates the exponentially scaled Airy function of the second kind defined in equation \ref{bie}.
80a3b2e0a   daniau   git-svn-id: https...
2184
2185
2186
2187
  %\begin{equation}
  % \label{bie}
  %bie(x)=Bi(x)\textrm{~if~}x\leq0 \qquad \qquad bie(x)=e^{-{2\over3}x^{3\over2}}Bi(x)\rm{~if~}x>0
  %\end{equation}
38581db0c   daniau   git-svn-id: https...
2188
2189
  \begin{equation}
   \label{bie}
80a3b2e0a   daniau   git-svn-id: https...
2190
2191
2192
2193
  bie(x)=\left\{\begin{array}{ll}
  Bi(x) & \textrm{~if~}x\leq0 \\
  e^{-{2\over3}x^{3\over2}}Bi(x) & \textrm{~if~}x>0
  \end{array}\right.
38581db0c   daniau   git-svn-id: https...
2194
2195
2196
2197
2198
2199
2200
2201
2202
2203
2204
2205
  \end{equation}
  
  Specific interfaces : \verb'bie,dbie'
  
  \subsubsection{aide}
  \begin{verbatim}
   aide(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates the exponentially scaled derivative of the Airy function as defined in equation \ref{aide}.
80a3b2e0a   daniau   git-svn-id: https...
2206
2207
2208
2209
  %\begin{equation}
  % \label{aide}
  %aie(x)=Ai^\prime(x) \textrm{~if~}x\leq0 \qquad \qquad aie(x)=e^{{2\over3}x^{3\over2}}Ai^\prime(x)\rm{~if~}x>0
  %\end{equation}
38581db0c   daniau   git-svn-id: https...
2210
2211
  \begin{equation}
   \label{aide}
80a3b2e0a   daniau   git-svn-id: https...
2212
2213
2214
2215
  aide(x)=\left\{\begin{array}{ll}
  Ai^\prime(x) & \textrm{~if~}x\leq0 \\
  e^{{2\over3}x^{3\over2}}Ai^\prime(x) & \textrm{~if~}x>0
  \end{array}\right.
38581db0c   daniau   git-svn-id: https...
2216
  \end{equation}
80a3b2e0a   daniau   git-svn-id: https...
2217

38581db0c   daniau   git-svn-id: https...
2218
2219
2220
2221
2222
2223
2224
2225
2226
2227
2228
2229
2230
  Specific interfaces : \verb'aide,daide'
  
  
  \subsubsection{bide}
  \begin{verbatim}
   bide(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates the exponentially scaled derivative of the Airy function of the second kind as defined in equation \ref{bide}.
  \begin{equation}
   \label{bide}
80a3b2e0a   daniau   git-svn-id: https...
2231
2232
2233
2234
  bie(x)=\left\{\begin{array}{ll}
  Bi^\prime(x) & \textrm{~if~}x\leq0 \\
  e^{-{2\over3}x^{3\over2}}Bi^\prime(x) & \textrm{~if~}x>0
  \end{array}\right.
38581db0c   daniau   git-svn-id: https...
2235
2236
2237
2238
2239
2240
2241
2242
2243
2244
2245
2246
2247
2248
2249
2250
2251
2252
2253
2254
2255
2256
2257
2258
2259
2260
2261
2262
2263
2264
2265
2266
2267
2268
2269
2270
2271
2272
2273
2274
2275
2276
2277
2278
2279
2280
2281
2282
  \end{equation}
  
  Specific interfaces : \verb'bide,dbide'
  
  
  \subsection{Miscellanous functions}
  \subsubsection{spenc}
  \begin{verbatim}
   spenc(x)
  \end{verbatim}
  \begin{itemize}
   \item x is a real
  \end{itemize}
  This function evaluates Spence function defined in equation \ref{spenc}.
  \begin{equation}
   \label{spenc}
  spenc(x)=s(x)=- \int_0^x {{ln(|1-t|)}\over t}dt
  \end{equation}
  
  Specific interfaces : \verb'spenc,dspenc'
  
  
  \subsubsection{inits}
  \begin{verbatim}
   inits(os,nos,eta)
  \end{verbatim}
  \begin{itemize}
   \item os is a real vector of length nos, containing the coefficients in an orthogonal series.
   \item nos is an integer
   \item eta is a real (Warning eta is a real(4) even with the double precision version) representing the requested accuracy.
  \end{itemize}
  This function initialize the orthogonal series so that inits is the number of terms needed to insure the error is no larger than eta.
  
  Specific interfaces : \verb'inits,initds'
  
  
  \subsubsection{csevl}
  \begin{verbatim}
   csevl(x,cs,n)
  \end{verbatim}
  \begin{itemize}
   \item x is a real in [-1,1]
   \item cs is a real vector of length n containing the coefficients of the Chebyshev serie.
   \item n is an integer
  \end{itemize}
  This function evaluates the Chebyshev series whose coefficients are stored in cs.
  
  Specific interfaces : \verb'csevl,dcsevl'
9158e74d6   daniau   git-svn-id: https...
2283
2284
2285
  
  
  \end{document}