legendre_rule.hpp 34 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486
  1. #ifndef _SCTL_LEGENDRE_RULE_HPP_
  2. #define _SCTL_LEGENDRE_RULE_HPP_
  3. # include <cstdlib>
  4. # include <cmath>
  5. # include <iostream>
  6. # include <fstream>
  7. # include <iomanip>
  8. # include <ctime>
  9. # include <cstring>
  10. # include <cassert>
  11. void cdgqf(int nt, int kind, double alpha, double beta, double t[], double wts[]);
  12. void cgqf(int nt, int kind, double alpha, double beta, double a, double b, double t[], double wts[]);
  13. double class_matrix(int kind, int m, double alpha, double beta, double aj[], double bj[]);
  14. void imtqlx(int n, double d[], double e[], double z[]);
  15. void parchk(int kind, int m, double alpha, double beta);
  16. double r8_abs(double x);
  17. double r8_epsilon();
  18. double r8_sign(double x);
  19. void r8mat_write(std::string output_filename, int m, int n, double table[]);
  20. void rule_write(int order, std::string filename, double x[], double w[], double r[]);
  21. 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);
  22. void sgqf(int nt, double aj[], double bj[], double zemu, double t[], double wts[]);
  23. //using namespace std;
  24. //****************************************************************************80
  25. /*
  26. int main ( int argc, char *argv[] )
  27. // ***************************************************************************80
  28. //
  29. // Purpose:
  30. //
  31. // MAIN is the main program for LEGENDRE_RULE.
  32. //
  33. // Discussion:
  34. //
  35. // This program computes a standard Gauss-Legendre quadrature rule
  36. // and writes it to a file.
  37. //
  38. // The user specifies:
  39. // * the ORDER (number of points) in the rule
  40. // * A, the left endpoint;
  41. // * B, the right endpoint;
  42. // * FILENAME, the root name of the output files.
  43. //
  44. // Licensing:
  45. //
  46. // This code is distributed under the GNU LGPL license.
  47. //
  48. // Modified:
  49. //
  50. // 23 February 2010
  51. //
  52. // Author:
  53. //
  54. // John Burkardt
  55. //
  56. {
  57. double a;
  58. double alpha;
  59. double b;
  60. double beta;
  61. string filename;
  62. int kind;
  63. int order;
  64. double *r;
  65. double *w;
  66. double *x;
  67. timestamp ( );
  68. std::cout << "\n";
  69. std::cout << "LEGENDRE_RULE\n";
  70. std::cout << " C++ version\n";
  71. std::cout << "\n";
  72. std::cout << " Compiled on " << __DATE__ << " at " << __TIME__ << ".\n";
  73. std::cout << "\n";
  74. std::cout << " Compute a Gauss-Legendre rule for approximating\n";
  75. std::cout << "\n";
  76. std::cout << " Integral ( A <= x <= B ) f(x) dx\n";
  77. std::cout << "\n";
  78. std::cout << " of order ORDER.\n";
  79. std::cout << "\n";
  80. std::cout << " The user specifies ORDER, A, B and FILENAME.\n";
  81. std::cout << "\n";
  82. std::cout << " Order is the number of points.\n";
  83. std::cout << "\n";
  84. std::cout << " A is the left endpoint.\n";
  85. std::cout << "\n";
  86. std::cout << " B is the right endpoint.\n";
  87. std::cout << "\n";
  88. std::cout << " FILENAME is used to generate 3 files:\n";
  89. std::cout << "\n";
  90. std::cout << " filename_w.txt - the weight file\n";
  91. std::cout << " filename_x.txt - the abscissa file.\n";
  92. std::cout << " filename_r.txt - the region file.\n";
  93. //
  94. // Initialize parameters;
  95. //
  96. alpha = 0.0;
  97. beta = 0.0;
  98. //
  99. // Get ORDER.
  100. //
  101. if ( 1 < argc )
  102. {
  103. order = atoi ( argv[1] );
  104. }
  105. else
  106. {
  107. std::cout << "\n";
  108. std::cout << " Enter the value of ORDER (1 or greater)\n";
  109. cin >> order;
  110. }
  111. //
  112. // Get A.
  113. //
  114. if ( 2 < argc )
  115. {
  116. a = atof ( argv[2] );
  117. }
  118. else
  119. {
  120. std::cout << "\n";
  121. std::cout << " Enter the value of A:\n";
  122. cin >> a;
  123. }
  124. //
  125. // Get B.
  126. //
  127. if ( 3 < argc )
  128. {
  129. b = atof ( argv[3] );
  130. }
  131. else
  132. {
  133. std::cout << "\n";
  134. std::cout << " Enter the value of B:\n";
  135. cin >> b;
  136. }
  137. //
  138. // Get FILENAME:
  139. //
  140. if ( 4 < argc )
  141. {
  142. filename = argv[4];
  143. }
  144. else
  145. {
  146. std::cout << "\n";
  147. std::cout << " Enter FILENAME, the \"root name\" of the quadrature files).\n";
  148. cin >> filename;
  149. }
  150. //
  151. // Input summary.
  152. //
  153. std::cout << "\n";
  154. std::cout << " ORDER = " << order << "\n";
  155. std::cout << " A = " << a << "\n";
  156. std::cout << " B = " << b << "\n";
  157. std::cout << " FILENAME = \"" << filename << "\".\n";
  158. //
  159. // Construct the rule.
  160. //
  161. w = new double[order];
  162. x = new double[order];
  163. kind = 1;
  164. cgqf ( order, kind, alpha, beta, a, b, x, w );
  165. //
  166. // Write the rule.
  167. //
  168. r = new double[2];
  169. r[0] = a;
  170. r[1] = b;
  171. rule_write ( order, filename, x, w, r );
  172. //
  173. // Terminate.
  174. //
  175. delete [] r;
  176. delete [] w;
  177. delete [] x;
  178. std::cout << "\n";
  179. std::cout << "LEGENDRE_RULE:\n";
  180. std::cout << " Normal end of execution.\n";
  181. std::cout << "\n";
  182. timestamp ( );
  183. return 0;
  184. }*/
  185. //****************************************************************************80
  186. inline void cdgqf ( int nt, int kind, double alpha, double beta, double t[],
  187. double wts[] )
  188. //****************************************************************************80
  189. //
  190. // Purpose:
  191. //
  192. // CDGQF computes a Gauss quadrature formula with default A, B and simple knots.
  193. //
  194. // Discussion:
  195. //
  196. // This routine computes all the knots and weights of a Gauss quadrature
  197. // formula with a classical weight function with default values for A and B,
  198. // and only simple knots.
  199. //
  200. // There are no moments checks and no printing is done.
  201. //
  202. // Use routine EIQFS to evaluate a quadrature computed by CGQFS.
  203. //
  204. // Licensing:
  205. //
  206. // This code is distributed under the GNU LGPL license.
  207. //
  208. // Modified:
  209. //
  210. // 08 January 2010
  211. //
  212. // Author:
  213. //
  214. // Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky.
  215. // C++ version by John Burkardt.
  216. //
  217. // Reference:
  218. //
  219. // Sylvan Elhay, Jaroslav Kautsky,
  220. // Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of
  221. // Interpolatory Quadrature,
  222. // ACM Transactions on Mathematical Software,
  223. // Volume 13, Number 4, December 1987, pages 399-415.
  224. //
  225. // Parameters:
  226. //
  227. // Input, int NT, the number of knots.
  228. //
  229. // Input, int KIND, the rule.
  230. // 1, Legendre, (a,b) 1.0
  231. // 2, Chebyshev Type 1, (a,b) ((b-x)*(x-a))^(-0.5)
  232. // 3, Gegenbauer, (a,b) ((b-x)*(x-a))^alpha
  233. // 4, Jacobi, (a,b) (b-x)^alpha*(x-a)^beta
  234. // 5, Generalized Laguerre, (a,inf) (x-a)^alpha*exp(-b*(x-a))
  235. // 6, Generalized Hermite, (-inf,inf) |x-a|^alpha*exp(-b*(x-a)^2)
  236. // 7, Exponential, (a,b) |x-(a+b)/2.0|^alpha
  237. // 8, Rational, (a,inf) (x-a)^alpha*(x+b)^beta
  238. //
  239. // Input, double ALPHA, the value of Alpha, if needed.
  240. //
  241. // Input, double BETA, the value of Beta, if needed.
  242. //
  243. // Output, double T[NT], the knots.
  244. //
  245. // Output, double WTS[NT], the weights.
  246. //
  247. {
  248. double *aj;
  249. double *bj;
  250. double zemu;
  251. parchk ( kind, 2 * nt, alpha, beta );
  252. //
  253. // Get the Jacobi matrix and zero-th moment.
  254. //
  255. aj = new double[nt];
  256. bj = new double[nt];
  257. zemu = class_matrix ( kind, nt, alpha, beta, aj, bj );
  258. //
  259. // Compute the knots and weights.
  260. //
  261. sgqf ( nt, aj, bj, zemu, t, wts );
  262. delete [] aj;
  263. delete [] bj;
  264. return;
  265. }
  266. //****************************************************************************80
  267. inline void cgqf ( int nt, int kind, double alpha, double beta, double a, double b,
  268. double t[], double wts[] )
  269. //****************************************************************************80
  270. //
  271. // Purpose:
  272. //
  273. // CGQF computes knots and weights of a Gauss quadrature formula.
  274. //
  275. // Discussion:
  276. //
  277. // The user may specify the interval (A,B).
  278. //
  279. // Only simple knots are produced.
  280. //
  281. // Use routine EIQFS to evaluate this quadrature formula.
  282. //
  283. // Licensing:
  284. //
  285. // This code is distributed under the GNU LGPL license.
  286. //
  287. // Modified:
  288. //
  289. // 16 February 2010
  290. //
  291. // Author:
  292. //
  293. // Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky.
  294. // C++ version by John Burkardt.
  295. //
  296. // Reference:
  297. //
  298. // Sylvan Elhay, Jaroslav Kautsky,
  299. // Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of
  300. // Interpolatory Quadrature,
  301. // ACM Transactions on Mathematical Software,
  302. // Volume 13, Number 4, December 1987, pages 399-415.
  303. //
  304. // Parameters:
  305. //
  306. // Input, int NT, the number of knots.
  307. //
  308. // Input, int KIND, the rule.
  309. // 1, Legendre, (a,b) 1.0
  310. // 2, Chebyshev Type 1, (a,b) ((b-x)*(x-a))^-0.5)
  311. // 3, Gegenbauer, (a,b) ((b-x)*(x-a))^alpha
  312. // 4, Jacobi, (a,b) (b-x)^alpha*(x-a)^beta
  313. // 5, Generalized Laguerre, (a,+oo) (x-a)^alpha*exp(-b*(x-a))
  314. // 6, Generalized Hermite, (-oo,+oo) |x-a|^alpha*exp(-b*(x-a)^2)
  315. // 7, Exponential, (a,b) |x-(a+b)/2.0|^alpha
  316. // 8, Rational, (a,+oo) (x-a)^alpha*(x+b)^beta
  317. // 9, Chebyshev Type 2, (a,b) ((b-x)*(x-a))^(+0.5)
  318. //
  319. // Input, double ALPHA, the value of Alpha, if needed.
  320. //
  321. // Input, double BETA, the value of Beta, if needed.
  322. //
  323. // Input, double A, B, the interval endpoints, or
  324. // other parameters.
  325. //
  326. // Output, double T[NT], the knots.
  327. //
  328. // Output, double WTS[NT], the weights.
  329. //
  330. {
  331. int i;
  332. int *mlt;
  333. int *ndx;
  334. //
  335. // Compute the Gauss quadrature formula for default values of A and B.
  336. //
  337. cdgqf ( nt, kind, alpha, beta, t, wts );
  338. //
  339. // Prepare to scale the quadrature formula to other weight function with
  340. // valid A and B.
  341. //
  342. mlt = new int[nt];
  343. for ( i = 0; i < nt; i++ )
  344. {
  345. mlt[i] = 1;
  346. }
  347. ndx = new int[nt];
  348. for ( i = 0; i < nt; i++ )
  349. {
  350. ndx[i] = i + 1;
  351. }
  352. scqf ( nt, t, mlt, wts, nt, ndx, wts, t, kind, alpha, beta, a, b );
  353. delete [] mlt;
  354. delete [] ndx;
  355. return;
  356. }
  357. //****************************************************************************80
  358. inline double class_matrix ( int kind, int m, double alpha, double beta, double aj[],
  359. double bj[] )
  360. //****************************************************************************80
  361. //
  362. // Purpose:
  363. //
  364. // CLASS_MATRIX computes the Jacobi matrix for a quadrature rule.
  365. //
  366. // Discussion:
  367. //
  368. // This routine computes the diagonal AJ and sub-diagonal BJ
  369. // elements of the order M tridiagonal symmetric Jacobi matrix
  370. // associated with the polynomials orthogonal with respect to
  371. // the weight function specified by KIND.
  372. //
  373. // For weight functions 1-7, M elements are defined in BJ even
  374. // though only M-1 are needed. For weight function 8, BJ(M) is
  375. // set to zero.
  376. //
  377. // The zero-th moment of the weight function is returned in ZEMU.
  378. //
  379. // Licensing:
  380. //
  381. // This code is distributed under the GNU LGPL license.
  382. //
  383. // Modified:
  384. //
  385. // 08 January 2010
  386. //
  387. // Author:
  388. //
  389. // Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky.
  390. // C++ version by John Burkardt.
  391. //
  392. // Reference:
  393. //
  394. // Sylvan Elhay, Jaroslav Kautsky,
  395. // Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of
  396. // Interpolatory Quadrature,
  397. // ACM Transactions on Mathematical Software,
  398. // Volume 13, Number 4, December 1987, pages 399-415.
  399. //
  400. // Parameters:
  401. //
  402. // Input, int KIND, the rule.
  403. // 1, Legendre, (a,b) 1.0
  404. // 2, Chebyshev Type 1, (a,b) ((b-x)*(x-a))^(-0.5)
  405. // 3, Gegenbauer, (a,b) ((b-x)*(x-a))^alpha
  406. // 4, Jacobi, (a,b) (b-x)^alpha*(x-a)^beta
  407. // 5, Generalized Laguerre, (a,inf) (x-a)^alpha*exp(-b*(x-a))
  408. // 6, Generalized Hermite, (-inf,inf) |x-a|^alpha*exp(-b*(x-a)^2)
  409. // 7, Exponential, (a,b) |x-(a+b)/2.0|^alpha
  410. // 8, Rational, (a,inf) (x-a)^alpha*(x+b)^beta
  411. //
  412. // Input, int M, the order of the Jacobi matrix.
  413. //
  414. // Input, double ALPHA, the value of Alpha, if needed.
  415. //
  416. // Input, double BETA, the value of Beta, if needed.
  417. //
  418. // Output, double AJ[M], BJ[M], the diagonal and subdiagonal
  419. // of the Jacobi matrix.
  420. //
  421. // Output, double CLASS_MATRIX, the zero-th moment.
  422. //
  423. {
  424. double a2b2;
  425. double ab;
  426. double aba;
  427. double abi;
  428. double abj;
  429. double abti;
  430. double apone;
  431. int i;
  432. double pi = 3.14159265358979323846264338327950;
  433. double zemu=0;
  434. parchk ( kind, 2 * m - 1, alpha, beta );
  435. /*
  436. double temp = r8_epsilon ( );
  437. double temp2 = 0.5;
  438. if ( 500.0 * temp < r8_abs ( pow ( tgamma ( temp2 ), 2 ) - pi ) )
  439. {
  440. std::cout << "\n";
  441. std::cout << "CLASS_MATRIX - Fatal error!\n";
  442. std::cout << " Gamma function does not match machine parameters.\n";
  443. exit ( 1 );
  444. }
  445. */
  446. if ( kind == 1 )
  447. {
  448. ab = 0.0;
  449. zemu = 2.0 / ( ab + 1.0 );
  450. for ( i = 0; i < m; i++ )
  451. {
  452. aj[i] = 0.0;
  453. }
  454. for ( i = 1; i <= m; i++ )
  455. {
  456. abi = i + ab * ( i % 2 );
  457. abj = 2 * i + ab;
  458. bj[i-1] = sqrt ( abi * abi / ( abj * abj - 1.0 ) );
  459. }
  460. }
  461. else if ( kind == 2 )
  462. {
  463. zemu = pi;
  464. for ( i = 0; i < m; i++ )
  465. {
  466. aj[i] = 0.0;
  467. }
  468. bj[0] = sqrt ( 0.5 );
  469. for ( i = 1; i < m; i++ )
  470. {
  471. bj[i] = 0.5;
  472. }
  473. }
  474. else if ( kind == 3 )
  475. {
  476. ab = alpha * 2.0;
  477. zemu = pow ( 2.0, ab + 1.0 ) * pow ( tgamma ( alpha + 1.0 ), 2 )
  478. / tgamma ( ab + 2.0 );
  479. for ( i = 0; i < m; i++ )
  480. {
  481. aj[i] = 0.0;
  482. }
  483. bj[0] = sqrt ( 1.0 / ( 2.0 * alpha + 3.0 ) );
  484. for ( i = 2; i <= m; i++ )
  485. {
  486. bj[i-1] = sqrt ( i * ( i + ab ) / ( 4.0 * pow ( i + alpha, 2 ) - 1.0 ) );
  487. }
  488. }
  489. else if ( kind == 4 )
  490. {
  491. ab = alpha + beta;
  492. abi = 2.0 + ab;
  493. zemu = pow ( 2.0, ab + 1.0 ) * tgamma ( alpha + 1.0 )
  494. * tgamma ( beta + 1.0 ) / tgamma ( abi );
  495. aj[0] = ( beta - alpha ) / abi;
  496. bj[0] = sqrt ( 4.0 * ( 1.0 + alpha ) * ( 1.0 + beta )
  497. / ( ( abi + 1.0 ) * abi * abi ) );
  498. a2b2 = beta * beta - alpha * alpha;
  499. for ( i = 2; i <= m; i++ )
  500. {
  501. abi = 2.0 * i + ab;
  502. aj[i-1] = a2b2 / ( ( abi - 2.0 ) * abi );
  503. abi = abi * abi;
  504. bj[i-1] = sqrt ( 4.0 * i * ( i + alpha ) * ( i + beta ) * ( i + ab )
  505. / ( ( abi - 1.0 ) * abi ) );
  506. }
  507. }
  508. else if ( kind == 5 )
  509. {
  510. zemu = tgamma ( alpha + 1.0 );
  511. for ( i = 1; i <= m; i++ )
  512. {
  513. aj[i-1] = 2.0 * i - 1.0 + alpha;
  514. bj[i-1] = sqrt ( i * ( i + alpha ) );
  515. }
  516. }
  517. else if ( kind == 6 )
  518. {
  519. zemu = tgamma ( ( alpha + 1.0 ) / 2.0 );
  520. for ( i = 0; i < m; i++ )
  521. {
  522. aj[i] = 0.0;
  523. }
  524. for ( i = 1; i <= m; i++ )
  525. {
  526. bj[i-1] = sqrt ( ( i + alpha * ( i % 2 ) ) / 2.0 );
  527. }
  528. }
  529. else if ( kind == 7 )
  530. {
  531. ab = alpha;
  532. zemu = 2.0 / ( ab + 1.0 );
  533. for ( i = 0; i < m; i++ )
  534. {
  535. aj[i] = 0.0;
  536. }
  537. for ( i = 1; i <= m; i++ )
  538. {
  539. abi = i + ab * ( i % 2 );
  540. abj = 2 * i + ab;
  541. bj[i-1] = sqrt ( abi * abi / ( abj * abj - 1.0 ) );
  542. }
  543. }
  544. else if ( kind == 8 )
  545. {
  546. ab = alpha + beta;
  547. zemu = tgamma ( alpha + 1.0 ) * tgamma ( - ( ab + 1.0 ) )
  548. / tgamma ( - beta );
  549. apone = alpha + 1.0;
  550. aba = ab * apone;
  551. aj[0] = - apone / ( ab + 2.0 );
  552. bj[0] = - aj[0] * ( beta + 1.0 ) / ( ab + 2.0 ) / ( ab + 3.0 );
  553. for ( i = 2; i <= m; i++ )
  554. {
  555. abti = ab + 2.0 * i;
  556. aj[i-1] = aba + 2.0 * ( ab + i ) * ( i - 1 );
  557. aj[i-1] = - aj[i-1] / abti / ( abti - 2.0 );
  558. }
  559. for ( i = 2; i <= m - 1; i++ )
  560. {
  561. abti = ab + 2.0 * i;
  562. bj[i-1] = i * ( alpha + i ) / ( abti - 1.0 ) * ( beta + i )
  563. / ( abti * abti ) * ( ab + i ) / ( abti + 1.0 );
  564. }
  565. bj[m-1] = 0.0;
  566. for ( i = 0; i < m; i++ )
  567. {
  568. bj[i] = sqrt ( bj[i] );
  569. }
  570. }
  571. return zemu;
  572. }
  573. //****************************************************************************80
  574. inline void imtqlx ( int n, double d[], double e[], double z[] )
  575. //****************************************************************************80
  576. //
  577. // Purpose:
  578. //
  579. // IMTQLX diagonalizes a symmetric tridiagonal matrix.
  580. //
  581. // Discussion:
  582. //
  583. // This routine is a slightly modified version of the EISPACK routine to
  584. // perform the implicit QL algorithm on a symmetric tridiagonal matrix.
  585. //
  586. // The authors thank the authors of EISPACK for permission to use this
  587. // routine.
  588. //
  589. // It has been modified to produce the product Q' * Z, where Z is an input
  590. // vector and Q is the orthogonal matrix diagonalizing the input matrix.
  591. // The changes consist (essentialy) of applying the orthogonal transformations
  592. // directly to Z as they are generated.
  593. //
  594. // Licensing:
  595. //
  596. // This code is distributed under the GNU LGPL license.
  597. //
  598. // Modified:
  599. //
  600. // 08 January 2010
  601. //
  602. // Author:
  603. //
  604. // Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky.
  605. // C++ version by John Burkardt.
  606. //
  607. // Reference:
  608. //
  609. // Sylvan Elhay, Jaroslav Kautsky,
  610. // Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of
  611. // Interpolatory Quadrature,
  612. // ACM Transactions on Mathematical Software,
  613. // Volume 13, Number 4, December 1987, pages 399-415.
  614. //
  615. // Roger Martin, James Wilkinson,
  616. // The Implicit QL Algorithm,
  617. // Numerische Mathematik,
  618. // Volume 12, Number 5, December 1968, pages 377-383.
  619. //
  620. // Parameters:
  621. //
  622. // Input, int N, the order of the matrix.
  623. //
  624. // Input/output, double D(N), the diagonal entries of the matrix.
  625. // On output, the information in D has been overwritten.
  626. //
  627. // Input/output, double E(N), the subdiagonal entries of the
  628. // matrix, in entries E(1) through E(N-1). On output, the information in
  629. // E has been overwritten.
  630. //
  631. // Input/output, double Z(N). On input, a vector. On output,
  632. // the value of Q' * Z, where Q is the matrix that diagonalizes the
  633. // input symmetric tridiagonal matrix.
  634. //
  635. {
  636. double b;
  637. double c;
  638. double f;
  639. double g;
  640. int i;
  641. int ii;
  642. int itn = 30;
  643. int j;
  644. int k;
  645. int l;
  646. int m=0;
  647. int mml;
  648. double p;
  649. double prec;
  650. double r;
  651. double s;
  652. prec = r8_epsilon ( );
  653. if ( n == 1 )
  654. {
  655. return;
  656. }
  657. e[n-1] = 0.0;
  658. for ( l = 1; l <= n; l++ )
  659. {
  660. j = 0;
  661. for ( ; ; )
  662. {
  663. for ( m = l; m <= n; m++ )
  664. {
  665. if ( m == n )
  666. {
  667. break;
  668. }
  669. if ( r8_abs ( e[m-1] ) <= prec * ( r8_abs ( d[m-1] ) + r8_abs ( d[m] ) ) )
  670. {
  671. break;
  672. }
  673. }
  674. p = d[l-1];
  675. if ( m == l )
  676. {
  677. break;
  678. }
  679. if ( itn <= j )
  680. {
  681. std::cout << "\n";
  682. std::cout << "IMTQLX - Fatal error!\n";
  683. std::cout << " Iteration limit exceeded\n";
  684. exit ( 1 );
  685. }
  686. j = j + 1;
  687. g = ( d[l] - p ) / ( 2.0 * e[l-1] );
  688. r = sqrt ( g * g + 1.0 );
  689. g = d[m-1] - p + e[l-1] / ( g + r8_abs ( r ) * r8_sign ( g ) );
  690. s = 1.0;
  691. c = 1.0;
  692. p = 0.0;
  693. mml = m - l;
  694. for ( ii = 1; ii <= mml; ii++ )
  695. {
  696. i = m - ii;
  697. f = s * e[i-1];
  698. b = c * e[i-1];
  699. if ( r8_abs ( g ) <= r8_abs ( f ) )
  700. {
  701. c = g / f;
  702. r = sqrt ( c * c + 1.0 );
  703. e[i] = f * r;
  704. s = 1.0 / r;
  705. c = c * s;
  706. }
  707. else
  708. {
  709. s = f / g;
  710. r = sqrt ( s * s + 1.0 );
  711. e[i] = g * r;
  712. c = 1.0 / r;
  713. s = s * c;
  714. }
  715. g = d[i] - p;
  716. r = ( d[i-1] - g ) * s + 2.0 * c * b;
  717. p = s * r;
  718. d[i] = g + p;
  719. g = c * r - b;
  720. f = z[i];
  721. z[i] = s * z[i-1] + c * f;
  722. z[i-1] = c * z[i-1] - s * f;
  723. }
  724. d[l-1] = d[l-1] - p;
  725. e[l-1] = g;
  726. e[m-1] = 0.0;
  727. }
  728. }
  729. //
  730. // Sorting.
  731. //
  732. for ( ii = 2; ii <= m; ii++ )
  733. {
  734. i = ii - 1;
  735. k = i;
  736. p = d[i-1];
  737. for ( j = ii; j <= n; j++ )
  738. {
  739. if ( d[j-1] < p )
  740. {
  741. k = j;
  742. p = d[j-1];
  743. }
  744. }
  745. if ( k != i )
  746. {
  747. d[k-1] = d[i-1];
  748. d[i-1] = p;
  749. p = z[i-1];
  750. z[i-1] = z[k-1];
  751. z[k-1] = p;
  752. }
  753. }
  754. return;
  755. }
  756. //****************************************************************************80
  757. inline void parchk ( int kind, int m, double alpha, double beta )
  758. //****************************************************************************80
  759. //
  760. // Purpose:
  761. //
  762. // PARCHK checks parameters ALPHA and BETA for classical weight functions.
  763. //
  764. // Licensing:
  765. //
  766. // This code is distributed under the GNU LGPL license.
  767. //
  768. // Modified:
  769. //
  770. // 07 January 2010
  771. //
  772. // Author:
  773. //
  774. // Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky.
  775. // C++ version by John Burkardt.
  776. //
  777. // Reference:
  778. //
  779. // Sylvan Elhay, Jaroslav Kautsky,
  780. // Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of
  781. // Interpolatory Quadrature,
  782. // ACM Transactions on Mathematical Software,
  783. // Volume 13, Number 4, December 1987, pages 399-415.
  784. //
  785. // Parameters:
  786. //
  787. // Input, int KIND, the rule.
  788. // 1, Legendre, (a,b) 1.0
  789. // 2, Chebyshev Type 1, (a,b) ((b-x)*(x-a))^(-0.5)
  790. // 3, Gegenbauer, (a,b) ((b-x)*(x-a))^alpha
  791. // 4, Jacobi, (a,b) (b-x)^alpha*(x-a)^beta
  792. // 5, Generalized Laguerre, (a,inf) (x-a)^alpha*exp(-b*(x-a))
  793. // 6, Generalized Hermite, (-inf,inf) |x-a|^alpha*exp(-b*(x-a)^2)
  794. // 7, Exponential, (a,b) |x-(a+b)/2.0|^alpha
  795. // 8, Rational, (a,inf) (x-a)^alpha*(x+b)^beta
  796. //
  797. // Input, int M, the order of the highest moment to
  798. // be calculated. This value is only needed when KIND = 8.
  799. //
  800. // Input, double ALPHA, BETA, the parameters, if required
  801. // by the value of KIND.
  802. //
  803. {
  804. double tmp;
  805. if ( kind <= 0 )
  806. {
  807. std::cout << "\n";
  808. std::cout << "PARCHK - Fatal error!\n";
  809. std::cout << " KIND <= 0.\n";
  810. exit ( 1 );
  811. }
  812. //
  813. // Check ALPHA for Gegenbauer, Jacobi, Laguerre, Hermite, Exponential.
  814. //
  815. if ( 3 <= kind && alpha <= -1.0 )
  816. {
  817. std::cout << "\n";
  818. std::cout << "PARCHK - Fatal error!\n";
  819. std::cout << " 3 <= KIND and ALPHA <= -1.\n";
  820. exit ( 1 );
  821. }
  822. //
  823. // Check BETA for Jacobi.
  824. //
  825. if ( kind == 4 && beta <= -1.0 )
  826. {
  827. std::cout << "\n";
  828. std::cout << "PARCHK - Fatal error!\n";
  829. std::cout << " KIND == 4 and BETA <= -1.0.\n";
  830. exit ( 1 );
  831. }
  832. //
  833. // Check ALPHA and BETA for rational.
  834. //
  835. if ( kind == 8 )
  836. {
  837. tmp = alpha + beta + m + 1.0;
  838. if ( 0.0 <= tmp || tmp <= beta )
  839. {
  840. std::cout << "\n";
  841. std::cout << "PARCHK - Fatal error!\n";
  842. std::cout << " KIND == 8 but condition on ALPHA and BETA fails.\n";
  843. exit ( 1 );
  844. }
  845. }
  846. return;
  847. }
  848. //****************************************************************************80
  849. inline double r8_abs ( double x )
  850. //****************************************************************************80
  851. //
  852. // Purpose:
  853. //
  854. // R8_ABS returns the absolute value of an R8.
  855. //
  856. // Licensing:
  857. //
  858. // This code is distributed under the GNU LGPL license.
  859. //
  860. // Modified:
  861. //
  862. // 14 November 2006
  863. //
  864. // Author:
  865. //
  866. // John Burkardt
  867. //
  868. // Parameters:
  869. //
  870. // Input, double X, the quantity whose absolute value is desired.
  871. //
  872. // Output, double R8_ABS, the absolute value of X.
  873. //
  874. {
  875. double value;
  876. if ( 0.0 <= x )
  877. {
  878. value = x;
  879. }
  880. else
  881. {
  882. value = - x;
  883. }
  884. return value;
  885. }
  886. //****************************************************************************80
  887. inline double r8_epsilon ( )
  888. //****************************************************************************80
  889. //
  890. // Purpose:
  891. //
  892. // R8_EPSILON returns the R8 roundoff unit.
  893. //
  894. // Discussion:
  895. //
  896. // The roundoff unit is a number R which is a power of 2 with the
  897. // property that, to the precision of the computer's arithmetic,
  898. // 1 < 1 + R
  899. // but
  900. // 1 = ( 1 + R / 2 )
  901. //
  902. // Licensing:
  903. //
  904. // This code is distributed under the GNU LGPL license.
  905. //
  906. // Modified:
  907. //
  908. // 01 July 2004
  909. //
  910. // Author:
  911. //
  912. // John Burkardt
  913. //
  914. // Parameters:
  915. //
  916. // Output, double R8_EPSILON, the R8 round-off unit.
  917. //
  918. {
  919. double value;
  920. value = 1.0;
  921. while ( 1.0 < ( double ) ( 1.0 + value ) )
  922. {
  923. value = value / 2.0;
  924. }
  925. value = 2.0 * value;
  926. return value;
  927. }
  928. //****************************************************************************80
  929. inline double r8_sign ( double x )
  930. //****************************************************************************80
  931. //
  932. // Purpose:
  933. //
  934. // R8_SIGN returns the sign of an R8.
  935. //
  936. // Licensing:
  937. //
  938. // This code is distributed under the GNU LGPL license.
  939. //
  940. // Modified:
  941. //
  942. // 18 October 2004
  943. //
  944. // Author:
  945. //
  946. // John Burkardt
  947. //
  948. // Parameters:
  949. //
  950. // Input, double X, the number whose sign is desired.
  951. //
  952. // Output, double R8_SIGN, the sign of X.
  953. //
  954. {
  955. double value;
  956. if ( x < 0.0 )
  957. {
  958. value = -1.0;
  959. }
  960. else
  961. {
  962. value = 1.0;
  963. }
  964. return value;
  965. }
  966. //****************************************************************************80
  967. inline void r8mat_write ( std::string output_filename, int m, int n, double table[] )
  968. //****************************************************************************80
  969. //
  970. // Purpose:
  971. //
  972. // R8MAT_WRITE writes an R8MAT file with no header.
  973. //
  974. // Licensing:
  975. //
  976. // This code is distributed under the GNU LGPL license.
  977. //
  978. // Modified:
  979. //
  980. // 29 June 2009
  981. //
  982. // Author:
  983. //
  984. // John Burkardt
  985. //
  986. // Parameters:
  987. //
  988. // Input, string OUTPUT_FILENAME, the output filename.
  989. //
  990. // Input, int M, the spatial dimension.
  991. //
  992. // Input, int N, the number of points.
  993. //
  994. // Input, double TABLE[M*N], the table data.
  995. //
  996. {
  997. int i;
  998. int j;
  999. std::ofstream output;
  1000. //
  1001. // Open the file.
  1002. //
  1003. output.open ( output_filename.c_str ( ) );
  1004. if ( !output )
  1005. {
  1006. std::cerr << "\n";
  1007. std::cerr << "R8MAT_WRITE - Fatal error!\n";
  1008. std::cerr << " Could not open the output file.\n";
  1009. return;
  1010. }
  1011. //
  1012. // Write the data.
  1013. //
  1014. for ( j = 0; j < n; j++ )
  1015. {
  1016. for ( i = 0; i < m; i++ )
  1017. {
  1018. output << " " << std::setw(24) << std::setprecision(16) << table[i+j*m];
  1019. }
  1020. output << "\n";
  1021. }
  1022. //
  1023. // Close the file.
  1024. //
  1025. output.close ( );
  1026. return;
  1027. }
  1028. //****************************************************************************80
  1029. inline void rule_write ( int order, std::string filename, double x[], double w[],
  1030. double r[] )
  1031. //****************************************************************************80
  1032. //
  1033. // Purpose:
  1034. //
  1035. // RULE_WRITE writes a quadrature rule to three files.
  1036. //
  1037. // Licensing:
  1038. //
  1039. // This code is distributed under the GNU LGPL license.
  1040. //
  1041. // Modified:
  1042. //
  1043. // 18 February 2010
  1044. //
  1045. // Author:
  1046. //
  1047. // John Burkardt
  1048. //
  1049. // Parameters:
  1050. //
  1051. // Input, int ORDER, the order of the rule.
  1052. //
  1053. // Input, double A, the left endpoint.
  1054. //
  1055. // Input, double B, the right endpoint.
  1056. //
  1057. // Input, string FILENAME, specifies the output filenames.
  1058. // "filename_w.txt", "filename_x.txt", "filename_r.txt"
  1059. // defining weights, abscissas, and region.
  1060. //
  1061. {
  1062. std::string filename_r;
  1063. std::string filename_w;
  1064. std::string filename_x;
  1065. filename_w = filename + "_w.txt";
  1066. filename_x = filename + "_x.txt";
  1067. filename_r = filename + "_r.txt";
  1068. std::cout << "\n";
  1069. std::cout << " Creating quadrature files.\n";
  1070. std::cout << "\n";
  1071. std::cout << " Root file name is \"" << filename << "\".\n";
  1072. std::cout << "\n";
  1073. std::cout << " Weight file will be \"" << filename_w << "\".\n";
  1074. std::cout << " Abscissa file will be \"" << filename_x << "\".\n";
  1075. std::cout << " Region file will be \"" << filename_r << "\".\n";
  1076. r8mat_write ( filename_w, 1, order, w );
  1077. r8mat_write ( filename_x, 1, order, x );
  1078. r8mat_write ( filename_r, 1, 2, r );
  1079. return;
  1080. }
  1081. //****************************************************************************80
  1082. inline void scqf ( int nt, double t[], int mlt[], double wts[], int nwts, int ndx[],
  1083. double swts[], double st[], int kind, double alpha, double beta, double a,
  1084. double b )
  1085. //****************************************************************************80
  1086. //
  1087. // Purpose:
  1088. //
  1089. // SCQF scales a quadrature formula to a nonstandard interval.
  1090. //
  1091. // Discussion:
  1092. //
  1093. // The arrays WTS and SWTS may coincide.
  1094. //
  1095. // The arrays T and ST may coincide.
  1096. //
  1097. // Licensing:
  1098. //
  1099. // This code is distributed under the GNU LGPL license.
  1100. //
  1101. // Modified:
  1102. //
  1103. // 16 February 2010
  1104. //
  1105. // Author:
  1106. //
  1107. // Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky.
  1108. // C++ version by John Burkardt.
  1109. //
  1110. // Reference:
  1111. //
  1112. // Sylvan Elhay, Jaroslav Kautsky,
  1113. // Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of
  1114. // Interpolatory Quadrature,
  1115. // ACM Transactions on Mathematical Software,
  1116. // Volume 13, Number 4, December 1987, pages 399-415.
  1117. //
  1118. // Parameters:
  1119. //
  1120. // Input, int NT, the number of knots.
  1121. //
  1122. // Input, double T[NT], the original knots.
  1123. //
  1124. // Input, int MLT[NT], the multiplicity of the knots.
  1125. //
  1126. // Input, double WTS[NWTS], the weights.
  1127. //
  1128. // Input, int NWTS, the number of weights.
  1129. //
  1130. // Input, int NDX[NT], used to index the array WTS.
  1131. // For more details see the comments in CAWIQ.
  1132. //
  1133. // Output, double SWTS[NWTS], the scaled weights.
  1134. //
  1135. // Output, double ST[NT], the scaled knots.
  1136. //
  1137. // Input, int KIND, the rule.
  1138. // 1, Legendre, (a,b) 1.0
  1139. // 2, Chebyshev Type 1, (a,b) ((b-x)*(x-a))^(-0.5)
  1140. // 3, Gegenbauer, (a,b) ((b-x)*(x-a))^alpha
  1141. // 4, Jacobi, (a,b) (b-x)^alpha*(x-a)^beta
  1142. // 5, Generalized Laguerre, (a,+oo) (x-a)^alpha*exp(-b*(x-a))
  1143. // 6, Generalized Hermite, (-oo,+oo) |x-a|^alpha*exp(-b*(x-a)^2)
  1144. // 7, Exponential, (a,b) |x-(a+b)/2.0|^alpha
  1145. // 8, Rational, (a,+oo) (x-a)^alpha*(x+b)^beta
  1146. // 9, Chebyshev Type 2, (a,b) ((b-x)*(x-a))^(+0.5)
  1147. //
  1148. // Input, double ALPHA, the value of Alpha, if needed.
  1149. //
  1150. // Input, double BETA, the value of Beta, if needed.
  1151. //
  1152. // Input, double A, B, the interval endpoints.
  1153. //
  1154. {
  1155. (void)(nwts); //SCTL_UNUSED(nwts);
  1156. double al=0;
  1157. double be=0;
  1158. int i;
  1159. int k;
  1160. int l;
  1161. double p;
  1162. double shft=0;
  1163. double slp=0;
  1164. double temp;
  1165. double tmp;
  1166. temp = r8_epsilon ( );
  1167. parchk ( kind, 1, alpha, beta );
  1168. if ( kind == 1 )
  1169. {
  1170. al = 0.0;
  1171. be = 0.0;
  1172. if ( r8_abs ( b - a ) <= temp )
  1173. {
  1174. std::cout << "\n";
  1175. std::cout << "SCQF - Fatal error!\n";
  1176. std::cout << " |B - A| too small.\n";
  1177. exit ( 1 );
  1178. }
  1179. shft = ( a + b ) / 2.0;
  1180. slp = ( b - a ) / 2.0;
  1181. }
  1182. else if ( kind == 2 )
  1183. {
  1184. al = -0.5;
  1185. be = -0.5;
  1186. if ( r8_abs ( b - a ) <= temp )
  1187. {
  1188. std::cout << "\n";
  1189. std::cout << "SCQF - Fatal error!\n";
  1190. std::cout << " |B - A| too small.\n";
  1191. exit ( 1 );
  1192. }
  1193. shft = ( a + b ) / 2.0;
  1194. slp = ( b - a ) / 2.0;
  1195. }
  1196. else if ( kind == 3 )
  1197. {
  1198. al = alpha;
  1199. be = alpha;
  1200. if ( r8_abs ( b - a ) <= temp )
  1201. {
  1202. std::cout << "\n";
  1203. std::cout << "SCQF - Fatal error!\n";
  1204. std::cout << " |B - A| too small.\n";
  1205. exit ( 1 );
  1206. }
  1207. shft = ( a + b ) / 2.0;
  1208. slp = ( b - a ) / 2.0;
  1209. }
  1210. else if ( kind == 4 )
  1211. {
  1212. al = alpha;
  1213. be = beta;
  1214. if ( r8_abs ( b - a ) <= temp )
  1215. {
  1216. std::cout << "\n";
  1217. std::cout << "SCQF - Fatal error!\n";
  1218. std::cout << " |B - A| too small.\n";
  1219. exit ( 1 );
  1220. }
  1221. shft = ( a + b ) / 2.0;
  1222. slp = ( b - a ) / 2.0;
  1223. }
  1224. else if ( kind == 5 )
  1225. {
  1226. if ( b <= 0.0 )
  1227. {
  1228. std::cout << "\n";
  1229. std::cout << "SCQF - Fatal error!\n";
  1230. std::cout << " B <= 0\n";
  1231. exit ( 1 );
  1232. }
  1233. shft = a;
  1234. slp = 1.0 / b;
  1235. al = alpha;
  1236. be = 0.0;
  1237. }
  1238. else if ( kind == 6 )
  1239. {
  1240. if ( b <= 0.0 )
  1241. {
  1242. std::cout << "\n";
  1243. std::cout << "SCQF - Fatal error!\n";
  1244. std::cout << " B <= 0.\n";
  1245. exit ( 1 );
  1246. }
  1247. shft = a;
  1248. slp = 1.0 / sqrt ( b );
  1249. al = alpha;
  1250. be = 0.0;
  1251. }
  1252. else if ( kind == 7 )
  1253. {
  1254. al = alpha;
  1255. be = 0.0;
  1256. if ( r8_abs ( b - a ) <= temp )
  1257. {
  1258. std::cout << "\n";
  1259. std::cout << "SCQF - Fatal error!\n";
  1260. std::cout << " |B - A| too small.\n";
  1261. exit ( 1 );
  1262. }
  1263. shft = ( a + b ) / 2.0;
  1264. slp = ( b - a ) / 2.0;
  1265. }
  1266. else if ( kind == 8 )
  1267. {
  1268. if ( a + b <= 0.0 )
  1269. {
  1270. std::cout << "\n";
  1271. std::cout << "SCQF - Fatal error!\n";
  1272. std::cout << " A + B <= 0.\n";
  1273. exit ( 1 );
  1274. }
  1275. shft = a;
  1276. slp = a + b;
  1277. al = alpha;
  1278. be = beta;
  1279. }
  1280. else if ( kind == 9 )
  1281. {
  1282. al = 0.5;
  1283. be = 0.5;
  1284. if ( r8_abs ( b - a ) <= temp )
  1285. {
  1286. std::cout << "\n";
  1287. std::cout << "SCQF - Fatal error!\n";
  1288. std::cout << " |B - A| too small.\n";
  1289. exit ( 1 );
  1290. }
  1291. shft = ( a + b ) / 2.0;
  1292. slp = ( b - a ) / 2.0;
  1293. }
  1294. p = pow ( slp, al + be + 1.0 );
  1295. for ( k = 0; k < nt; k++ )
  1296. {
  1297. st[k] = shft + slp * t[k];
  1298. l = abs ( ndx[k] );
  1299. if ( l != 0 )
  1300. {
  1301. tmp = p;
  1302. for ( i = l - 1; i <= l - 1 + mlt[k] - 1; i++ )
  1303. {
  1304. swts[i] = wts[i] * tmp;
  1305. tmp = tmp * slp;
  1306. }
  1307. }
  1308. }
  1309. return;
  1310. }
  1311. //****************************************************************************80
  1312. inline void sgqf ( int nt, double aj[], double bj[], double zemu, double t[],
  1313. double wts[] )
  1314. //****************************************************************************80
  1315. //
  1316. // Purpose:
  1317. //
  1318. // SGQF computes knots and weights of a Gauss Quadrature formula.
  1319. //
  1320. // Discussion:
  1321. //
  1322. // This routine computes all the knots and weights of a Gauss quadrature
  1323. // formula with simple knots from the Jacobi matrix and the zero-th
  1324. // moment of the weight function, using the Golub-Welsch technique.
  1325. //
  1326. // Licensing:
  1327. //
  1328. // This code is distributed under the GNU LGPL license.
  1329. //
  1330. // Modified:
  1331. //
  1332. // 08 January 2010
  1333. //
  1334. // Author:
  1335. //
  1336. // Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky.
  1337. // C++ version by John Burkardt.
  1338. //
  1339. // Reference:
  1340. //
  1341. // Sylvan Elhay, Jaroslav Kautsky,
  1342. // Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of
  1343. // Interpolatory Quadrature,
  1344. // ACM Transactions on Mathematical Software,
  1345. // Volume 13, Number 4, December 1987, pages 399-415.
  1346. //
  1347. // Parameters:
  1348. //
  1349. // Input, int NT, the number of knots.
  1350. //
  1351. // Input, double AJ[NT], the diagonal of the Jacobi matrix.
  1352. //
  1353. // Input/output, double BJ[NT], the subdiagonal of the Jacobi
  1354. // matrix, in entries 1 through NT-1. On output, BJ has been overwritten.
  1355. //
  1356. // Input, double ZEMU, the zero-th moment of the weight function.
  1357. //
  1358. // Output, double T[NT], the knots.
  1359. //
  1360. // Output, double WTS[NT], the weights.
  1361. //
  1362. {
  1363. int i;
  1364. //
  1365. // Exit if the zero-th moment is not positive.
  1366. //
  1367. if ( zemu <= 0.0 )
  1368. {
  1369. std::cout << "\n";
  1370. std::cout << "SGQF - Fatal error!\n";
  1371. std::cout << " ZEMU <= 0.\n";
  1372. exit ( 1 );
  1373. }
  1374. //
  1375. // Set up vectors for IMTQLX.
  1376. //
  1377. for ( i = 0; i < nt; i++ )
  1378. {
  1379. t[i] = aj[i];
  1380. }
  1381. wts[0] = sqrt ( zemu );
  1382. for ( i = 1; i < nt; i++ )
  1383. {
  1384. wts[i] = 0.0;
  1385. }
  1386. //
  1387. // Diagonalize the Jacobi matrix.
  1388. //
  1389. imtqlx ( nt, t, bj, wts );
  1390. for ( i = 0; i < nt; i++ )
  1391. {
  1392. wts[i] = wts[i] * wts[i];
  1393. }
  1394. return;
  1395. }
  1396. //****************************************************************************80
  1397. #endif