legendre_rule.hpp 34 KB

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