1 #ifndef _LEGENDRE_RULE_HPP_ 2 #define _LEGENDRE_RULE_HPP_ 6 void cdgqf(
int nt,
int kind,
double alpha,
double beta,
double t[],
double wts[]);
7 void cgqf(
int nt,
int kind,
double alpha,
double beta,
double a,
double b,
double t[],
double wts[]);
8 double class_matrix(
int kind,
int m,
double alpha,
double beta,
double aj[],
double bj[]);
9 void imtqlx(
int n,
double d[],
double e[],
double z[]);
10 void parchk(
int kind,
int m,
double alpha,
double beta);
11 double r8_abs(
double x);
13 double r8_sign(
double x);
14 void r8mat_write(std::string output_filename,
int m,
int n,
double table[]);
15 void rule_write(
int order, std::string filename,
double x[],
double w[],
double r[]);
16 void scqf(
int nt,
double t[],
int mlt[],
double wts[],
int nwts,
int ndx[],
double swts[],
double st[],
int kind,
double alpha,
double beta,
double a,
double b);
17 void sgqf(
int nt,
double aj[],
double bj[],
double zemu,
double t[],
double wts[]);
28 # include <pvfmm/legendre_rule.hpp> 201 inline void cdgqf (
int nt,
int kind,
double alpha,
double beta,
double t[],
268 parchk ( kind, 2 * nt, alpha, beta );
275 zemu = class_matrix ( kind, nt, alpha, beta, aj, bj );
279 sgqf ( nt, aj, bj, zemu, t, wts );
288 inline void cgqf (
int nt,
int kind,
double alpha,
double beta,
double a,
double b,
289 double t[],
double wts[] )
359 cdgqf ( nt, kind, alpha, beta, t, wts );
365 for ( i = 0; i < nt; i++ )
370 for ( i = 0; i < nt; i++ )
374 scqf ( nt, t, mlt, wts, nt, ndx, wts, t, kind, alpha, beta, a, b );
383 inline double class_matrix (
int kind,
int m,
double alpha,
double beta,
double aj[],
458 double pi = 3.14159265358979323846264338327950;
461 parchk ( kind, 2 * m - 1, alpha, beta );
478 zemu = 2.0 / ( ab + 1.0 );
480 for ( i = 0; i < m; i++ )
485 for ( i = 1; i <= m; i++ )
487 abi = i + ab * ( i % 2 );
489 bj[i-1] = sqrt ( abi * abi / ( abj * abj - 1.0 ) );
492 else if ( kind == 2 )
496 for ( i = 0; i < m; i++ )
501 bj[0] = sqrt ( 0.5 );
502 for ( i = 1; i < m; i++ )
507 else if ( kind == 3 )
510 zemu = pow ( 2.0, ab + 1.0 ) * pow ( tgamma ( alpha + 1.0 ), 2 )
511 / tgamma ( ab + 2.0 );
513 for ( i = 0; i < m; i++ )
518 bj[0] = sqrt ( 1.0 / ( 2.0 * alpha + 3.0 ) );
519 for ( i = 2; i <= m; i++ )
521 bj[i-1] = sqrt ( i * ( i + ab ) / ( 4.0 * pow ( i + alpha, 2 ) - 1.0 ) );
524 else if ( kind == 4 )
528 zemu = pow ( 2.0, ab + 1.0 ) * tgamma ( alpha + 1.0 )
529 * tgamma ( beta + 1.0 ) / tgamma ( abi );
530 aj[0] = ( beta - alpha ) / abi;
531 bj[0] = sqrt ( 4.0 * ( 1.0 + alpha ) * ( 1.0 + beta )
532 / ( ( abi + 1.0 ) * abi * abi ) );
533 a2b2 = beta * beta - alpha * alpha;
535 for ( i = 2; i <= m; i++ )
538 aj[i-1] = a2b2 / ( ( abi - 2.0 ) * abi );
540 bj[i-1] = sqrt ( 4.0 * i * ( i + alpha ) * ( i + beta ) * ( i + ab )
541 / ( ( abi - 1.0 ) * abi ) );
544 else if ( kind == 5 )
546 zemu = tgamma ( alpha + 1.0 );
548 for ( i = 1; i <= m; i++ )
550 aj[i-1] = 2.0 * i - 1.0 + alpha;
551 bj[i-1] = sqrt ( i * ( i + alpha ) );
554 else if ( kind == 6 )
556 zemu = tgamma ( ( alpha + 1.0 ) / 2.0 );
558 for ( i = 0; i < m; i++ )
563 for ( i = 1; i <= m; i++ )
565 bj[i-1] = sqrt ( ( i + alpha * ( i % 2 ) ) / 2.0 );
568 else if ( kind == 7 )
571 zemu = 2.0 / ( ab + 1.0 );
573 for ( i = 0; i < m; i++ )
578 for ( i = 1; i <= m; i++ )
580 abi = i + ab * ( i % 2 );
582 bj[i-1] = sqrt ( abi * abi / ( abj * abj - 1.0 ) );
585 else if ( kind == 8 )
588 zemu = tgamma ( alpha + 1.0 ) * tgamma ( - ( ab + 1.0 ) )
592 aj[0] = - apone / ( ab + 2.0 );
593 bj[0] = - aj[0] * ( beta + 1.0 ) / ( ab + 2.0 ) / ( ab + 3.0 );
594 for ( i = 2; i <= m; i++ )
597 aj[i-1] = aba + 2.0 * ( ab + i ) * ( i - 1 );
598 aj[i-1] = - aj[i-1] / abti / ( abti - 2.0 );
601 for ( i = 2; i <= m - 1; i++ )
604 bj[i-1] = i * ( alpha + i ) / ( abti - 1.0 ) * ( beta + i )
605 / ( abti * abti ) * ( ab + i ) / ( abti + 1.0 );
608 for ( i = 0; i < m; i++ )
610 bj[i] = sqrt ( bj[i] );
618 inline void imtqlx (
int n,
double d[],
double e[],
double z[] )
698 prec = r8_epsilon ( );
707 for ( l = 1; l <= n; l++ )
712 for ( m = l; m <= n; m++ )
719 if ( r8_abs ( e[m-1] ) <= prec * ( r8_abs ( d[m-1] ) + r8_abs ( d[m] ) ) )
732 std::cout <<
"IMTQLX - Fatal error!\n";
733 std::cout <<
" Iteration limit exceeded\n";
737 g = ( d[l] - p ) / ( 2.0 * e[l-1] );
738 r = sqrt ( g * g + 1.0 );
739 g = d[m-1] - p + e[l-1] / ( g + r8_abs ( r ) * r8_sign ( g ) );
745 for ( ii = 1; ii <= mml; ii++ )
751 if ( r8_abs ( g ) <= r8_abs ( f ) )
754 r = sqrt ( c * c + 1.0 );
762 r = sqrt ( s * s + 1.0 );
768 r = ( d[i-1] - g ) * s + 2.0 * c * b;
773 z[i] = s * z[i-1] + c * f;
774 z[i-1] = c * z[i-1] - s * f;
784 for ( ii = 2; ii <= m; ii++ )
790 for ( j = ii; j <= n; j++ )
812 inline void parchk (
int kind,
int m,
double alpha,
double beta )
865 std::cout <<
"PARCHK - Fatal error!\n";
866 std::cout <<
" KIND <= 0.\n";
872 if ( 3 <= kind && alpha <= -1.0 )
875 std::cout <<
"PARCHK - Fatal error!\n";
876 std::cout <<
" 3 <= KIND and ALPHA <= -1.\n";
882 if ( kind == 4 && beta <= -1.0 )
885 std::cout <<
"PARCHK - Fatal error!\n";
886 std::cout <<
" KIND == 4 and BETA <= -1.0.\n";
894 tmp = alpha + beta + m + 1.0;
895 if ( 0.0 <= tmp || tmp <= beta )
898 std::cout <<
"PARCHK - Fatal error!\n";
899 std::cout <<
" KIND == 8 but condition on ALPHA and BETA fails.\n";
907 inline double r8_abs (
double x )
948 inline double r8_epsilon ( )
985 while ( 1.0 < (
double ) ( 1.0 + value ) )
996 inline double r8_sign (
double x )
1037 inline void r8mat_write ( std::string output_filename,
int m,
int n,
double table[] )
1070 std::ofstream output;
1074 output.open ( output_filename.c_str ( ) );
1079 std::cerr <<
"R8MAT_WRITE - Fatal error!\n";
1080 std::cerr <<
" Could not open the output file.\n";
1086 for ( j = 0; j < n; j++ )
1088 for ( i = 0; i < m; i++ )
1090 output <<
" " << std::setw(24) << std::setprecision(16) << table[i+j*m];
1103 inline void rule_write (
int order, std::string filename,
double x[],
double w[],
1137 std::string filename_r;
1138 std::string filename_w;
1139 std::string filename_x;
1141 filename_w = filename +
"_w.txt";
1142 filename_x = filename +
"_x.txt";
1143 filename_r = filename +
"_r.txt";
1146 std::cout <<
" Creating quadrature files.\n";
1148 std::cout <<
" Root file name is \"" << filename <<
"\".\n";
1150 std::cout <<
" Weight file will be \"" << filename_w <<
"\".\n";
1151 std::cout <<
" Abscissa file will be \"" << filename_x <<
"\".\n";
1152 std::cout <<
" Region file will be \"" << filename_r <<
"\".\n";
1154 r8mat_write ( filename_w, 1, order, w );
1155 r8mat_write ( filename_x, 1, order, x );
1156 r8mat_write ( filename_r, 1, 2, r );
1162 inline void scqf (
int nt,
double t[],
int mlt[],
double wts[],
int nwts,
int ndx[],
1163 double swts[],
double st[],
int kind,
double alpha,
double beta,
double a,
1248 temp = r8_epsilon ( );
1250 parchk ( kind, 1, alpha, beta );
1256 if ( r8_abs ( b - a ) <= temp )
1259 std::cout <<
"SCQF - Fatal error!\n";
1260 std::cout <<
" |B - A| too small.\n";
1263 shft = ( a + b ) / 2.0;
1264 slp = ( b - a ) / 2.0;
1266 else if ( kind == 2 )
1270 if ( r8_abs ( b - a ) <= temp )
1273 std::cout <<
"SCQF - Fatal error!\n";
1274 std::cout <<
" |B - A| too small.\n";
1277 shft = ( a + b ) / 2.0;
1278 slp = ( b - a ) / 2.0;
1280 else if ( kind == 3 )
1284 if ( r8_abs ( b - a ) <= temp )
1287 std::cout <<
"SCQF - Fatal error!\n";
1288 std::cout <<
" |B - A| too small.\n";
1291 shft = ( a + b ) / 2.0;
1292 slp = ( b - a ) / 2.0;
1294 else if ( kind == 4 )
1299 if ( r8_abs ( b - a ) <= temp )
1302 std::cout <<
"SCQF - Fatal error!\n";
1303 std::cout <<
" |B - A| too small.\n";
1306 shft = ( a + b ) / 2.0;
1307 slp = ( b - a ) / 2.0;
1309 else if ( kind == 5 )
1314 std::cout <<
"SCQF - Fatal error!\n";
1315 std::cout <<
" B <= 0\n";
1323 else if ( kind == 6 )
1328 std::cout <<
"SCQF - Fatal error!\n";
1329 std::cout <<
" B <= 0.\n";
1333 slp = 1.0 / sqrt ( b );
1337 else if ( kind == 7 )
1341 if ( r8_abs ( b - a ) <= temp )
1344 std::cout <<
"SCQF - Fatal error!\n";
1345 std::cout <<
" |B - A| too small.\n";
1348 shft = ( a + b ) / 2.0;
1349 slp = ( b - a ) / 2.0;
1351 else if ( kind == 8 )
1356 std::cout <<
"SCQF - Fatal error!\n";
1357 std::cout <<
" A + B <= 0.\n";
1365 else if ( kind == 9 )
1369 if ( r8_abs ( b - a ) <= temp )
1372 std::cout <<
"SCQF - Fatal error!\n";
1373 std::cout <<
" |B - A| too small.\n";
1376 shft = ( a + b ) / 2.0;
1377 slp = ( b - a ) / 2.0;
1380 p = pow ( slp, al + be + 1.0 );
1382 for ( k = 0; k < nt; k++ )
1384 st[k] = shft + slp * t[k];
1390 for ( i = l - 1; i <= l - 1 + mlt[k] - 1; i++ )
1392 swts[i] = wts[i] * tmp;
1401 inline void sgqf (
int nt,
double aj[],
double bj[],
double zemu,
double t[],
1460 std::cout <<
"SGQF - Fatal error!\n";
1461 std::cout <<
" ZEMU <= 0.\n";
1467 for ( i = 0; i < nt; i++ )
1471 wts[0] = sqrt ( zemu );
1472 for ( i = 1; i < nt; i++ )
1479 imtqlx ( nt, t, bj, wts );
1481 for ( i = 0; i < nt; i++ )
1483 wts[i] = wts[i] * wts[i];
1490 inline void timestamp ( )
1519 # define TIME_SIZE 40 1521 static char time_buffer[TIME_SIZE];
1522 const struct std::tm *tm_ptr;
1526 now = std::time ( NULL );
1527 tm_ptr = std::localtime ( &now );
1529 len = std::strftime ( time_buffer, TIME_SIZE,
"%d %B %Y %I:%M:%S %p", tm_ptr );
1531 std::cout << time_buffer <<
"\n";