Cloudy
Spectral Synthesis Code for Astrophysics
Loading...
Searching...
No Matches
thirdparty.h
Go to the documentation of this file.
1/* This file contains routines (perhaps in modified form) by third parties.
2 * Use and distribution of these works are determined by their respective copyrights. */
3
4#ifndef THIRDPARTY_H_
5#define THIRDPARTY_H_
6
7#include "physconst.h"
8
9/*============================================================================*/
10
11/* these are the routines in thirdparty.cpp */
12
13bool linfit(
14 long n,
15 const double x[], /* x[n] */
16 const double y[], /* y[n] */
17 double &a,
18 double &siga,
19 double &b,
20 double &sigb
21);
22
25static const int NPRE_FACTORIAL = 171;
26
27/* largest value of fct function corresponding to factorial in six j sjs calculation */
28static const int MXDSF = 340;
29
31double factorial(long n);
32
34double lfactorial(long n);
35
36complex<double> cdgamma(complex<double> x);
37
38double bessel_j0(double x);
39double bessel_y0(double x);
40double bessel_j1(double x);
41double bessel_y1(double x);
42double bessel_jn(int n, double x);
43double bessel_yn(int n, double x);
44
45double ellpk(double x);
46
51double expn(int n, double x);
52
54double erfce(double);
55
56double igam (double a, double x);
57double igamc (double a, double x);
58double igamc_scaled(double a, double x);
59
60double bessel_k0(double x);
61double bessel_k0_scaled(double x);
62double bessel_k1(double x);
63double bessel_k1_scaled(double x);
64void bessel_k0_k1(double x, double* k0val, double* k1val);
65void bessel_k0_k1_scaled(double x, double* k0val, double* k1val);
66
67double bessel_i0(double x);
68double bessel_i0_scaled(double x);
69double bessel_i1(double x);
70double bessel_i1_scaled(double x);
71void bessel_i0_i1(double x, double* k0val, double* k1val);
72void bessel_i0_i1_scaled(double x, double* k0val, double* k1val);
73
74#ifndef HAVE_SINCOS
75// this is a GNU extension
76inline void sincos(double x, double* s, double* c)
77{
78 *s = sin(x);
79 *c = cos(x);
80}
81#endif
82
85double e1(double x);
86
88double e1_scaled(double x);
89
92double e2(double t);
93
94/* random number generator written by David Blackman and Sebastiano Vigna (vigna@acm.org) */
95void xoroshiro128plus(uint64* pool, size_t size, uint64 state[], size_t ns);
96/* this call is equivalent to skipping ahead 2^64 random variates in the sequence */
97void xoroshiro128plus_jump(uint64& state0, uint64& state1);
98/* another random number generator written by David Blackman and Sebastiano Vigna (vigna@acm.org) */
99void xoshiro256starstar(uint64* pool, size_t size, uint64 state[], size_t ns);
100/* this call is equivalent to skipping ahead 2^128 random variates in the sequence */
101void xoshiro256starstar_jump(uint64& state0, uint64& state1, uint64& state2, uint64& state3);
102
103/* simple random number generator written by Sebastiano Vigna (vigna@acm.org) */
104uint64 splitmix64(uint64& state);
105
106/*============================================================================*/
107
108/* these are the routines in thirdparty_interpolate.cpp */
109
110void spline_cubic_set( long n, const double t[], const double y[], double ypp[],
111 int ibcbeg, double ybcbeg, int ibcend, double ybcend );
112void spline_cubic_val( long n, const double t[], double tval, const double y[], const double ypp[],
113 double *yval, double *ypval, double *yppval );
114
127inline void spline(const double x[],
128 const double y[],
129 long int n,
130 double yp1,
131 double ypn,
132 double y2a[])
133{
134 int ibcbeg = ( yp1 > 0.99999e30 ) ? 2 : 1;
135 double ybcbeg = ( yp1 > 0.99999e30 ) ? 0. : yp1;
136 int ibcend = ( ypn > 0.99999e30 ) ? 2 : 1;
137 double ybcend = ( ypn > 0.99999e30 ) ? 0. : ypn;
138 spline_cubic_set( n, x, y, y2a, ibcbeg, ybcbeg, ibcend, ybcend );
139 return;
140}
141
150inline void splint(const double xa[],
151 const double ya[],
152 const double y2a[],
153 long int n,
154 double x,
155 double *y)
156{
157 spline_cubic_val( n, xa, x, ya, y2a, y, NULL, NULL );
158 return;
159}
160
161inline void spldrv(const double xa[],
162 const double ya[],
163 const double y2a[],
164 long int n,
165 double x,
166 double *y)
167{
168 spline_cubic_val( n, xa, x, ya, y2a, NULL, y, NULL );
169 return;
170}
171
172/* wrapper routine for splint that checks whether x-value is within bounds
173 * if the x-value is out of bounds, a flag will be raised and the function
174 * will be evaluated at the nearest boundary */
175/* >>chng 03 jan 15, added splint_safe, PvH */
176inline void splint_safe(const double xa[],
177 const double ya[],
178 const double y2a[],
179 long int n,
180 double x,
181 double *y,
182 bool *lgOutOfBounds)
183{
184 double xsafe;
185
186 const double lo_bound = MIN2(xa[0],xa[n-1]);
187 const double hi_bound = MAX2(xa[0],xa[n-1]);
188 const double SAFETY = MAX2(hi_bound-lo_bound,1.)*10.*DBL_EPSILON;
189
190 DEBUG_ENTRY( "splint_safe()" );
191
192 if( x < lo_bound-SAFETY )
193 {
194 xsafe = lo_bound;
195 *lgOutOfBounds = true;
196 }
197 else if( x > hi_bound+SAFETY )
198 {
199 xsafe = hi_bound;
200 *lgOutOfBounds = true;
201 }
202 else
203 {
204 xsafe = x;
205 *lgOutOfBounds = false;
206 }
207
208 splint(xa,ya,y2a,n,xsafe,y);
209 return;
210}
211
212/* wrapper routine for spldrv that checks whether x-value is within bounds
213 * if the x-value is out of bounds, a flag will be raised and the function
214 * will be evaluated at the nearest boundary */
215/* >>chng 03 jan 15, added spldrv_safe, PvH */
216inline void spldrv_safe(const double xa[],
217 const double ya[],
218 const double y2a[],
219 long int n,
220 double x,
221 double *y,
222 bool *lgOutOfBounds)
223{
224 double xsafe;
225
226 const double lo_bound = MIN2(xa[0],xa[n-1]);
227 const double hi_bound = MAX2(xa[0],xa[n-1]);
228 const double SAFETY = MAX2(fabs(lo_bound),fabs(hi_bound))*10.*DBL_EPSILON;
229
230 DEBUG_ENTRY( "spldrv_safe()" );
231
232 if( x < lo_bound-SAFETY )
233 {
234 xsafe = lo_bound;
235 *lgOutOfBounds = true;
236 }
237 else if( x > hi_bound+SAFETY )
238 {
239 xsafe = hi_bound;
240 *lgOutOfBounds = true;
241 }
242 else
243 {
244 xsafe = x;
245 *lgOutOfBounds = false;
246 }
247
248 spldrv(xa,ya,y2a,n,xsafe,y);
249 return;
250}
251
260double lagrange(const double x[], /* x[n] */
261 const double y[], /* y[n] */
262 long n,
263 double xval);
264
267template<class T>
268inline long hunt_bisect(const T x[], /* x[n] */
269 long n,
270 T xval)
271{
272 /* do bisection hunt */
273 long ilo = 0, ihi = n-1;
274 while( ihi-ilo > 1 )
275 {
276 long imid = (ilo+ihi)/2;
277 if( xval < x[imid] )
278 ihi = imid;
279 else
280 ilo = imid;
281 }
282 return ilo;
283}
284
286template<class T>
287inline long hunt_bisect(const vector<T>& x,
288 T xval)
289{
290 return hunt_bisect(get_ptr(x), x.size(), xval);
291}
292
293template<class T>
294inline long hunt_bisect(const vector<T>& x,
295 size_t n,
296 T xval)
297{
298 ASSERT( n <= x.size() );
299 return hunt_bisect(get_ptr(x), n, xval);
300}
301
304template<class T>
305inline long hunt_bisect_reverse(const T x[], /* x[n] */
306 long n,
307 T xval)
308{
309 /* do bisection hunt */
310 long ilo = 0, ihi = n-1;
311 while( ihi-ilo > 1 )
312 {
313 long imid = (ilo+ihi)/2;
314 if( xval <= x[imid] )
315 ilo = imid;
316 else
317 ihi = imid;
318 }
319 return ilo;
320}
321
323template<class T>
324inline long hunt_bisect_reverse(const vector<T>& x,
325 T xval)
326{
327 return hunt_bisect_reverse(get_ptr(x), x.size(), xval);
328}
329
330template<class T>
331inline long hunt_bisect_reverse(const vector<T>& x,
332 size_t n,
333 T xval)
334{
335 ASSERT( n <= x.size() );
336 return hunt_bisect_reverse(get_ptr(x), n, xval);
337}
338
346template<class T>
347T linint(const T x[], /* x[n] */
348 const T y[], /* y[n] */
349 long n,
350 T xval)
351{
352 T yval;
353
354 ASSERT( n >= 2 );
355
356 if( xval <= x[0] )
357 yval = y[0];
358 else if( xval >= x[n-1] )
359 yval = y[n-1];
360 else
361 {
362 long ilo = hunt_bisect( x, n, xval );
363 T deriv = (y[ilo+1]-y[ilo])/(x[ilo+1]-x[ilo]);
364 yval = y[ilo] + deriv*(xval-x[ilo]);
365 }
366 return yval;
367}
368
369/*============================================================================*/
370
371/* these are the routines in thirdparty_lapack.cpp */
372
373/* there are wrappers for lapack linear algebra routines.
374 * there are two versions of the lapack routines - a fortran
375 * version that I converted to C with forc to use if nothing else is available
376 * (included in the Cloudy distribution),
377 * and an option to link into an external lapack library that may be optimized
378 * for your machine. By default the tralated C routines will be used.
379 * To use your machine's lapack library instead, define the macro
380 * LAPACK and link into your library. This is usually done with a command line
381 * switch "-DLAPACK" on the compiler command, and the linker option "-llapack"
382 */
383
392void getrf_wrapper(long M, long N, double *A, long lda, int32 *ipiv, int32 *info);
393
405void getrs_wrapper(char trans, long N, long nrhs, double *A, long lda, int32 *ipiv, double *B, long ldb, int32 *info);
406
407double dawson(double x);
408
409void humlik(int n, const realnum x[], realnum y, realnum k[]);
410
412void FastVoigtH(realnum a, const realnum v[], realnum y[], size_t n);
413
414// calculates y[i] = H(a,v[i]) as defined in Eq. 9-44 of Mihalas
415inline void VoigtH(realnum a, const realnum v[], realnum y[], int n)
416{
417 if( a <= 0.1f )
418 {
419 FastVoigtH( a, v, y, n );
420 }
421 else
422 {
423 humlik( n, v, a, y );
424 }
425}
426
427// calculates y[i] = U(a,v[i]) as defined in Eq. 9-45 of Mihalas
428inline void VoigtU(realnum a, const realnum v[], realnum y[], int n)
429{
430 VoigtH( a, v, y, n );
431 for( int i=0; i < n; ++i )
432 y[i] /= realnum(SQRTPI);
433}
434
435// VoigtH0(a) returns the value H(a,0) following Eq. 9-44 of Mihalas
436inline double VoigtH0(double a)
437{
438 return erfce(a);
439}
440
441// VoigtU0(a) returns the value U(a,0) following Eq. 9-45 of Mihalas
442inline double VoigtU0(double a)
443{
444 return erfce(a)/SQRTPI;
445}
446
448static const unsigned int NMD5 = 32;
449
451string MD5file(const char* fnam, access_scheme scheme=AS_DEFAULT);
453string MD5datafile(const char* fnam, access_scheme scheme=AS_DEFAULT);
455string MD5datastream(fstream& ioFile);
457string MD5string(const string& str);
458void MD5string(const string& str, uint64 md5sum[2]);
459
460void test_expn();
461void chbfit(double a, double b, vector<double>& c, double (*func)(double));
462
463void svd(const int nm, const int m, const int n, double *a,
464 double *w, bool matu, double *u, bool matv,
465 double *v, int *ierr, double *rv1);
466
467// Evaluate the Gegenbauer (aka ultraspherical) polynomial C_n^(alpha)(x)
468double gegenbauer(long n, double al, double x);
469
470// Diagonal generating function for ultraspherical polynomials,
471// following Badnell analysis. Integer argument to constructor is
472// sum of raised and lowered index to the ultraspherical function.
474{
475 long m_a, m_n;
476 double m_x, m_c1, m_c2;
477public:
478UltraGen(long a, double x) : m_a(a), m_n(0), m_x(x), m_c1(0.), m_c2(1.)
479 {}
480 void step()
481 {
482 double c1 = m_c1;
483 --m_a;
484 m_c1 = 2.*m_a*(m_c2-m_x*c1)/(m_n+2*m_a);
485 m_c2 = 2.*m_a*(m_x*m_c2-c1)/(m_n+1);
486 ++m_n;
487 }
488 double val() const
489 {
490 return m_c2;
491 }
492};
493
494// this is the triangular inequality for angular momenta in 2*J format
495inline bool Triangle2(long J1, long J2, long J3)
496{
497 return ( J1 >= 0 && J2 >= 0 && J3 >= 0 &&
498 J1 >= abs(J2-J3) && J1 <= J2+J3 &&
499 J1%2 == (J2+J3)%2 );
500}
501
502// 6j Wigner evaluation, original routine Fortran Nigel Badnell 6j for Autostructure
503double sjs(long j1, long j2, long j3, long l1, long l2, long l3);
504// version of the above with less range in j, is exported for testing purposes
505double SixJFull(long j1, long j2, long j3, long j4, long j5, long j6);
506// Implementation utlizing recursion relation for vector of sixj values,
507// should be robust to high j
508void rec6j(double *sixcof, double l2, double l3,
509 double l4, double l5, double l6, double *l1min,
510 double *l1max, double *lmatch, long ndim, long *ier);
511
512// the following code was taken (and altered) from:
513// http://www.geeksforgeeks.org/find-all-shortest-unique-prefixes-to-represent-each-word-in-a-given-list/
514static const size_t TRIESZ = 128;
515
517{
519 size_t freq; // To store frequency
521 {
522 memset(child, 0, TRIESZ*sizeof(trieNode*));
523 freq = 0;
524 }
525 trieNode(const trieNode&) = delete;
526 trieNode& operator= (const trieNode&) = delete;
528 {
529 for( size_t i=0; i < TRIESZ; i++ )
530 delete child[i];
531 }
532};
533
534void insertToken(trieNode* root, const string& token);
535size_t findUniqueLen(trieNode* root, const string& token);
536
537// the following code was taken (and altered) from:
538// https://en.wikipedia.org/wiki/Levenshtein_distance
539
540size_t LevenshteinDistance(const string& s, const string& t);
541
542#endif /* THIRDPARTY_H_ */
#define NULL
Definition cddefines.h:115
#define ASSERT(exp)
Definition cddefines.h:637
float realnum
Definition cddefines.h:127
#define MAX2(a, b)
Definition cddefines.h:828
T * get_ptr(T *v)
Definition cddefines.h:1115
#define MIN2(a, b)
Definition cddefines.h:807
#define DEBUG_ENTRY(funcname)
Definition cddefines.h:730
double m_x
Definition thirdparty.h:476
long m_a
Definition thirdparty.h:475
void step()
Definition thirdparty.h:480
long m_n
Definition thirdparty.h:475
double m_c1
Definition thirdparty.h:476
double val() const
Definition thirdparty.h:488
double m_c2
Definition thirdparty.h:476
UltraGen(long a, double x)
Definition thirdparty.h:478
access_scheme
Definition cpu.h:257
@ AS_DEFAULT
Definition cpu.h:257
static const double SAFETY
Definition grains_qheat.cpp:54
Definition thirdparty.h:517
trieNode & operator=(const trieNode &)=delete
trieNode * child[TRIESZ]
Definition thirdparty.h:518
trieNode()
Definition thirdparty.h:520
trieNode(const trieNode &)=delete
~trieNode()
Definition thirdparty.h:527
size_t freq
Definition thirdparty.h:519
void xoshiro256starstar(uint64 *pool, size_t size, uint64 state[], size_t ns)
Definition thirdparty.cpp:3779
void xoroshiro128plus(uint64 *pool, size_t size, uint64 state[], size_t ns)
Definition thirdparty.cpp:3712
double igamc_scaled(double a, double x)
Definition thirdparty.cpp:2101
void sincos(double x, double *s, double *c)
Definition thirdparty.h:76
double erfce(double)
Definition thirdparty.cpp:1975
double bessel_j1(double x)
Definition thirdparty.cpp:1297
double igam(double a, double x)
Definition thirdparty.cpp:2128
double bessel_k1_scaled(double x)
Definition thirdparty.cpp:2541
double factorial(long n)
Definition thirdparty.cpp:357
double bessel_k0(double x)
Definition thirdparty.cpp:2472
double ellpk(double x)
Definition thirdparty.cpp:1677
double dawson(double x)
Definition thirdparty.cpp:3988
double bessel_i1_scaled(double x)
Definition thirdparty.cpp:2727
double bessel_yn(int n, double x)
Definition thirdparty.cpp:1535
void spldrv(const double xa[], const double ya[], const double y2a[], long int n, double x, double *y)
Definition thirdparty.h:161
double bessel_y0(double x)
Definition thirdparty.cpp:1104
double bessel_k0_scaled(double x)
Definition thirdparty.cpp:2495
static const unsigned int NMD5
Definition thirdparty.h:448
size_t findUniqueLen(trieNode *root, const string &token)
Definition thirdparty.cpp:5153
void humlik(int n, const realnum x[], realnum y, realnum k[])
Definition thirdparty.cpp:4163
void getrs_wrapper(char trans, long N, long nrhs, double *A, long lda, int32 *ipiv, double *B, long ldb, int32 *info)
Definition thirdparty_lapack.cpp:69
size_t LevenshteinDistance(const string &s, const string &t)
Definition thirdparty.cpp:5173
double igamc(double a, double x)
Definition thirdparty.cpp:2083
double bessel_jn(int n, double x)
Definition thirdparty.cpp:1400
double VoigtH0(double a)
Definition thirdparty.h:436
bool Triangle2(long J1, long J2, long J3)
Definition thirdparty.h:495
long hunt_bisect_reverse(const T x[], long n, T xval)
Definition thirdparty.h:305
string MD5string(const string &str)
Definition thirdparty.cpp:4437
void rec6j(double *sixcof, double l2, double l3, double l4, double l5, double l6, double *l1min, double *l1max, double *lmatch, long ndim, long *ier)
Definition thirdparty.cpp:3380
double lfactorial(long n)
Definition thirdparty.cpp:757
double lagrange(const double x[], const double y[], long n, double xval)
Definition thirdparty_interpolate.cpp:430
string MD5datastream(fstream &ioFile)
Definition thirdparty.cpp:4360
double gegenbauer(long n, double al, double x)
Definition thirdparty.cpp:3236
double e1_scaled(double x)
Definition thirdparty.cpp:3068
void spline_cubic_set(long n, const double t[], const double y[], double ypp[], int ibcbeg, double ybcbeg, int ibcend, double ybcend)
Definition thirdparty_interpolate.cpp:105
void spldrv_safe(const double xa[], const double ya[], const double y2a[], long int n, double x, double *y, bool *lgOutOfBounds)
Definition thirdparty.h:216
void test_expn()
Definition thirdparty.cpp:3210
complex< double > cdgamma(complex< double > x)
Definition thirdparty.cpp:790
void bessel_i0_i1(double x, double *k0val, double *k1val)
Definition thirdparty.cpp:2748
void VoigtH(realnum a, const realnum v[], realnum y[], int n)
Definition thirdparty.h:415
long hunt_bisect(const T x[], long n, T xval)
Definition thirdparty.h:268
static const int NPRE_FACTORIAL
Definition thirdparty.h:25
void splint_safe(const double xa[], const double ya[], const double y2a[], long int n, double x, double *y, bool *lgOutOfBounds)
Definition thirdparty.h:176
void spline_cubic_val(long n, const double t[], double tval, const double y[], const double ypp[], double *yval, double *ypval, double *yppval)
Definition thirdparty_interpolate.cpp:328
void svd(const int nm, const int m, const int n, double *a, double *w, bool matu, double *u, bool matv, double *v, int *ierr, double *rv1)
Definition thirdparty.cpp:4669
double SixJFull(long j1, long j2, long j3, long j4, long j5, long j6)
Definition thirdparty.cpp:3330
double bessel_i1(double x)
Definition thirdparty.cpp:2706
double bessel_i0_scaled(double x)
Definition thirdparty.cpp:2685
realnum FastVoigtH(realnum a, realnum v)
Definition thirdparty.cpp:4027
void bessel_k0_k1(double x, double *k0val, double *k1val)
Definition thirdparty.cpp:2564
void bessel_k0_k1_scaled(double x, double *k0val, double *k1val)
Definition thirdparty.cpp:2592
double bessel_j0(double x)
Definition thirdparty.cpp:1066
double VoigtU0(double a)
Definition thirdparty.h:442
double bessel_i0(double x)
Definition thirdparty.cpp:2664
double e1(double x)
Definition thirdparty.cpp:3030
T linint(const T x[], const T y[], long n, T xval)
Definition thirdparty.h:347
static const size_t TRIESZ
Definition thirdparty.h:514
double e2(double t)
Definition thirdparty.cpp:3100
void bessel_i0_i1_scaled(double x, double *k0val, double *k1val)
Definition thirdparty.cpp:2775
void spline(const double x[], const double y[], long int n, double yp1, double ypn, double y2a[])
Definition thirdparty.h:127
void chbfit(double a, double b, vector< double > &c, double(*func)(double))
Definition thirdparty.cpp:3189
double bessel_k1(double x)
Definition thirdparty.cpp:2518
void getrf_wrapper(long M, long N, double *A, long lda, int32 *ipiv, int32 *info)
Definition thirdparty_lapack.cpp:47
string MD5file(const char *fnam, access_scheme scheme=AS_DEFAULT)
Definition thirdparty.cpp:4331
void insertToken(trieNode *root, const string &token)
Definition thirdparty.cpp:5129
void splint(const double xa[], const double ya[], const double y2a[], long int n, double x, double *y)
Definition thirdparty.h:150
double expn(int n, double x)
Definition thirdparty.cpp:1757
void xoroshiro128plus_jump(uint64 &state0, uint64 &state1)
Definition thirdparty.cpp:3738
string MD5datafile(const char *fnam, access_scheme scheme=AS_DEFAULT)
Definition thirdparty.cpp:4350
void xoshiro256starstar_jump(uint64 &state0, uint64 &state1, uint64 &state2, uint64 &state3)
Definition thirdparty.cpp:3811
bool linfit(long n, const double x[], const double y[], double &a, double &siga, double &b, double &sigb)
Definition thirdparty.cpp:47
double sjs(long j1, long j2, long j3, long l1, long l2, long l3)
Definition thirdparty.cpp:3266
double bessel_y1(double x)
Definition thirdparty.cpp:1324
static const int MXDSF
Definition thirdparty.h:28
void VoigtU(realnum a, const realnum v[], realnum y[], int n)
Definition thirdparty.h:428
uint64 splitmix64(uint64 &state)
Definition thirdparty.cpp:3870