Intrepid
Intrepid_PolylibDef.hpp
Go to the documentation of this file.
1 /*
2 // @HEADER
3 // ************************************************************************
4 //
5 // Intrepid Package
6 // Copyright (2007) Sandia Corporation
7 //
8 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
9 // license for use of this work by or on behalf of the U.S. Government.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Pavel Bochev (pbboche@sandia.gov)
39 // Denis Ridzal (dridzal@sandia.gov), or
40 // Kara Peterson (kjpeter@sandia.gov)
41 //
42 // ************************************************************************
43 // @HEADER
44 */
45 
47 //
48 // File: Intrepid_PolylibDef.hpp
49 //
50 // For more information, please see: http://www.nektar.info
51 //
52 // The MIT License
53 //
54 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
55 // Department of Aeronautics, Imperial College London (UK), and Scientific
56 // Computing and Imaging Institute, University of Utah (USA).
57 //
58 // License for the specific language governing rights and limitations under
59 // Permission is hereby granted, free of charge, to any person obtaining a
60 // copy of this software and associated documentation files (the "Software"),
61 // to deal in the Software without restriction, including without limitation
62 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
63 // and/or sell copies of the Software, and to permit persons to whom the
64 // Software is furnished to do so, subject to the following conditions:
65 //
66 // The above copyright notice and this permission notice shall be included
67 // in all copies or substantial portions of the Software.
68 //
69 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
70 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
71 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
72 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
73 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
74 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
75 // DEALINGS IN THE SOFTWARE.
76 //
77 // Description:
78 // This file is redistributed with the Intrepid package. It should be used
79 // in accordance with the above MIT license, at the request of the authors.
80 // This file is NOT covered by the usual Intrepid/Trilinos LGPL license.
81 //
82 // Origin: Nektar++ library, http://www.nektar.info, downloaded on
83 // March 10, 2009.
84 //
86 
87 
94 #ifdef _MSC_VER
95 #include "winmath.h"
96 #endif
97 
98 
99 namespace Intrepid {
100 
102 #define INTREPID_POLYLIB_STOP 50
103 
105 #define INTREPID_POLYLIB_POLYNOMIAL_DEFLATION 0
106 
107 #ifdef INTREPID_POLYLIB_POLYNOMIAL_DEFLATION
108 #define jacobz(n,z,alpha,beta) Jacobz(n,z,alpha,beta)
110 #else
111 #define jacobz(n,z,alpha,beta) JacZeros(n,z,alpha,beta)
113 #endif
114 
115 
116 template <class Scalar>
117 void IntrepidPolylib::zwgj (Scalar *z, Scalar *w, const int np, const Scalar alpha, const Scalar beta){
118  int i;
119  Scalar fac, one = 1.0, two = 2.0, apb = alpha + beta;
120 
121  IntrepidPolylib::jacobz (np,z,alpha,beta);
122  IntrepidPolylib::jacobd (np,z,w,np,alpha,beta);
123 
124  fac = std::pow(two,apb + one)*IntrepidPolylib::gammaF(alpha + np + one)*IntrepidPolylib::gammaF(beta + np + one);
125  fac /= IntrepidPolylib::gammaF((Scalar)(np + one))*IntrepidPolylib::gammaF(apb + np + one);
126 
127  for(i = 0; i < np; ++i) w[i] = fac/(w[i]*w[i]*(one-z[i]*z[i]));
128 
129  return;
130 }
131 
132 
133 template <class Scalar>
134 void IntrepidPolylib::zwgrjm(Scalar *z, Scalar *w, const int np, const Scalar alpha, const Scalar beta){
135 
136  if(np == 1){
137  z[0] = 0.0;
138  w[0] = 2.0;
139  }
140  else{
141  int i;
142  Scalar fac, one = 1.0, two = 2.0, apb = alpha + beta;
143 
144  z[0] = -one;
145  IntrepidPolylib::jacobz (np-1,z+1,alpha,beta+1);
146  IntrepidPolylib::jacobfd (np,z,w,(Scalar*)0,np-1,alpha,beta);
147 
148  fac = std::pow(two,apb)*IntrepidPolylib::gammaF(alpha + np)*IntrepidPolylib::gammaF(beta + np);
149  fac /= IntrepidPolylib::gammaF((Scalar)np)*(beta + np)*IntrepidPolylib::gammaF(apb + np + 1);
150 
151  for(i = 0; i < np; ++i) w[i] = fac*(1-z[i])/(w[i]*w[i]);
152  w[0] *= (beta + one);
153  }
154 
155  return;
156 }
157 
158 
159 template <class Scalar>
160 void IntrepidPolylib::zwgrjp(Scalar *z, Scalar *w, const int np, const Scalar alpha, const Scalar beta){
161 
162  if(np == 1){
163  z[0] = 0.0;
164  w[0] = 2.0;
165  }
166  else{
167  int i;
168  Scalar fac, one = 1.0, two = 2.0, apb = alpha + beta;
169 
170  IntrepidPolylib::jacobz (np-1,z,alpha+1,beta);
171  z[np-1] = one;
172  IntrepidPolylib::jacobfd (np,z,w,(Scalar*)0,np-1,alpha,beta);
173 
174  fac = std::pow(two,apb)*IntrepidPolylib::gammaF(alpha + np)*IntrepidPolylib::gammaF(beta + np);
175  fac /= IntrepidPolylib::gammaF((Scalar)np)*(alpha + np)*IntrepidPolylib::gammaF(apb + np + 1);
176 
177  for(i = 0; i < np; ++i) w[i] = fac*(1+z[i])/(w[i]*w[i]);
178  w[np-1] *= (alpha + one);
179  }
180 
181  return;
182 }
183 
184 
185 template <class Scalar>
186 void IntrepidPolylib::zwglj(Scalar *z, Scalar *w, const int np, const Scalar alpha, const Scalar beta){
187 
188  if( np == 1 ){
189  z[0] = 0.0;
190  w[0] = 2.0;
191  }
192  else{
193  int i;
194  Scalar fac, one = 1.0, apb = alpha + beta, two = 2.0;
195 
196  z[0] = -one;
197  z[np-1] = one;
198  IntrepidPolylib::jacobz (np-2,z + 1,alpha + one,beta + one);
199  IntrepidPolylib::jacobfd (np,z,w,(Scalar*)0,np-1,alpha,beta);
200 
201  fac = std::pow(two,apb + 1)*IntrepidPolylib::gammaF(alpha + np)*IntrepidPolylib::gammaF(beta + np);
202  fac /= (np-1)*IntrepidPolylib::gammaF((Scalar)np)*IntrepidPolylib::gammaF(alpha + beta + np + one);
203 
204  for(i = 0; i < np; ++i) w[i] = fac/(w[i]*w[i]);
205  w[0] *= (beta + one);
206  w[np-1] *= (alpha + one);
207  }
208 
209  return;
210 }
211 
212 
213 template <class Scalar>
214 void IntrepidPolylib::Dgj(Scalar *D, const Scalar *z, const int np, const Scalar alpha, const Scalar beta)
215 {
216 
217  Scalar one = 1.0, two = 2.0;
218 
219  if (np <= 0){
220  D[0] = 0.0;
221  }
222  else{
223  int i,j;
224  Scalar *pd;
225 
226  pd = (Scalar *)malloc(np*sizeof(Scalar));
227  IntrepidPolylib::jacobd(np,z,pd,np,alpha,beta);
228 
229  for (i = 0; i < np; i++){
230  for (j = 0; j < np; j++){
231 
232  if (i != j)
233  //D[i*np+j] = pd[j]/(pd[i]*(z[j]-z[i])); <--- This is either a bug, or the derivative matrix is not defined consistently.
234  D[i*np+j] = pd[i]/(pd[j]*(z[i]-z[j]));
235  else
236  D[i*np+j] = (alpha - beta + (alpha + beta + two)*z[j])/
237  (two*(one - z[j]*z[j]));
238  }
239  }
240  free(pd);
241  }
242  return;
243 }
244 
245 
246 template <class Scalar>
247 void IntrepidPolylib::Dgrjm(Scalar *D, const Scalar *z, const int np, const Scalar alpha, const Scalar beta)
248 {
249 
250  if (np <= 0){
251  D[0] = 0.0;
252  }
253  else{
254  int i, j;
255  Scalar one = 1.0, two = 2.0;
256  Scalar *pd;
257 
258  pd = (Scalar *)malloc(np*sizeof(Scalar));
259 
260  pd[0] = std::pow(-one,np-1)*IntrepidPolylib::gammaF((Scalar)np+beta+one);
261  pd[0] /= IntrepidPolylib::gammaF((Scalar)np)*IntrepidPolylib::gammaF(beta+two);
262  IntrepidPolylib::jacobd(np-1,z+1,pd+1,np-1,alpha,beta+1);
263  for(i = 1; i < np; ++i) pd[i] *= (1+z[i]);
264 
265  for (i = 0; i < np; i++) {
266  for (j = 0; j < np; j++){
267  if (i != j)
268  //D[i*np+j] = pd[j]/(pd[i]*(z[j]-z[i])); <--- This is either a bug, or the derivative matrix is not defined consistently.
269  D[i*np+j] = pd[i]/(pd[j]*(z[i]-z[j]));
270  else {
271  if(j == 0)
272  D[i*np+j] = -(np + alpha + beta + one)*(np - one)/
273  (two*(beta + two));
274  else
275  D[i*np+j] = (alpha - beta + one + (alpha + beta + one)*z[j])/
276  (two*(one - z[j]*z[j]));
277  }
278  }
279  }
280  free(pd);
281  }
282  return;
283 }
284 
285 
286 template <class Scalar>
287 void IntrepidPolylib::Dgrjp(Scalar *D, const Scalar *z, const int np, const Scalar alpha, const Scalar beta)
288 {
289 
290  if (np <= 0){
291  D[0] = 0.0;
292  }
293  else{
294  int i, j;
295  Scalar one = 1.0, two = 2.0;
296  Scalar *pd;
297 
298  pd = (Scalar *)malloc(np*sizeof(Scalar));
299 
300 
301  IntrepidPolylib::jacobd(np-1,z,pd,np-1,alpha+1,beta);
302  for(i = 0; i < np-1; ++i) pd[i] *= (1-z[i]);
303  pd[np-1] = -IntrepidPolylib::gammaF((Scalar)np+alpha+one);
304  pd[np-1] /= IntrepidPolylib::gammaF((Scalar)np)*IntrepidPolylib::gammaF(alpha+two);
305 
306  for (i = 0; i < np; i++) {
307  for (j = 0; j < np; j++){
308  if (i != j)
309  //D[i*np+j] = pd[j]/(pd[i]*(z[j]-z[i])); <--- This is either a bug, or the derivative matrix is not defined consistently.
310  D[i*np+j] = pd[i]/(pd[j]*(z[i]-z[j]));
311  else {
312  if(j == np-1)
313  D[i*np+j] = (np + alpha + beta + one)*(np - one)/
314  (two*(alpha + two));
315  else
316  D[i*np+j] = (alpha - beta - one + (alpha + beta + one)*z[j])/
317  (two*(one - z[j]*z[j]));
318  }
319  }
320  }
321  free(pd);
322  }
323  return;
324 }
325 
326 
327 template <class Scalar>
328 void IntrepidPolylib::Dglj(Scalar *D, const Scalar *z, const int np, const Scalar alpha, const Scalar beta)
329 {
330 
331  if (np <= 0){
332  D[0] = 0.0;
333  }
334  else{
335  int i, j;
336  Scalar one = 1.0, two = 2.0;
337  Scalar *pd;
338 
339  pd = (Scalar *)malloc(np*sizeof(Scalar));
340 
341  pd[0] = two*std::pow(-one,np)*IntrepidPolylib::gammaF((Scalar)np + beta);
342  pd[0] /= IntrepidPolylib::gammaF((Scalar)np - one)*IntrepidPolylib::gammaF(beta + two);
343  IntrepidPolylib::jacobd(np-2,z+1,pd+1,np-2,alpha+1,beta+1);
344  for(i = 1; i < np-1; ++i) pd[i] *= (one-z[i]*z[i]);
345  pd[np-1] = -two*IntrepidPolylib::gammaF((Scalar)np + alpha);
346  pd[np-1] /= IntrepidPolylib::gammaF((Scalar)np - one)*IntrepidPolylib::gammaF(alpha + two);
347 
348  for (i = 0; i < np; i++) {
349  for (j = 0; j < np; j++){
350  if (i != j)
351  //D[i*np+j] = pd[j]/(pd[i]*(z[j]-z[i])); <--- This is either a bug, or the derivative matrix is not defined consistently.
352  D[i*np+j] = pd[i]/(pd[j]*(z[i]-z[j]));
353  else {
354  if (j == 0)
355  D[i*np+j] = (alpha - (np-1)*(np + alpha + beta))/(two*(beta+ two));
356  else if (j == np-1)
357  D[i*np+j] =-(beta - (np-1)*(np + alpha + beta))/(two*(alpha+ two));
358  else
359  D[i*np+j] = (alpha - beta + (alpha + beta)*z[j])/
360  (two*(one - z[j]*z[j]));
361  }
362  }
363  }
364  free(pd);
365  }
366  return;
367 }
368 
369 
370 template <class Scalar>
371 Scalar IntrepidPolylib::hgj (const int i, const Scalar z, const Scalar *zgj,
372  const int np, const Scalar alpha, const Scalar beta)
373 {
374 
375  Scalar zi, dz, p, pd, h;
376 
377  zi = *(zgj+i);
378  dz = z - zi;
379  if (std::abs(dz) < INTREPID_TOL) return 1.0;
380 
381  IntrepidPolylib::jacobd (1, &zi, &pd , np, alpha, beta);
382  IntrepidPolylib::jacobfd(1, &z , &p, (Scalar*)0 , np, alpha, beta);
383  h = p/(pd*dz);
384 
385  return h;
386 }
387 
388 
389 template <class Scalar>
390 Scalar IntrepidPolylib::hgrjm (const int i, const Scalar z, const Scalar *zgrj,
391  const int np, const Scalar alpha, const Scalar beta)
392 {
393 
394  Scalar zi, dz, p, pd, h;
395 
396  zi = *(zgrj+i);
397  dz = z - zi;
398  if (std::abs(dz) < INTREPID_TOL) return 1.0;
399 
400  IntrepidPolylib::jacobfd (1, &zi, &p , (Scalar*)0, np-1, alpha, beta + 1);
401  // need to use this routine in case zi = -1 or 1
402  IntrepidPolylib::jacobd (1, &zi, &pd, np-1, alpha, beta + 1);
403  h = (1.0 + zi)*pd + p;
404  IntrepidPolylib::jacobfd (1, &z, &p, (Scalar*)0, np-1, alpha, beta + 1);
405  h = (1.0 + z )*p/(h*dz);
406 
407  return h;
408 }
409 
410 
411 template <class Scalar>
412 Scalar IntrepidPolylib::hgrjp (const int i, const Scalar z, const Scalar *zgrj,
413  const int np, const Scalar alpha, const Scalar beta)
414 {
415 
416  Scalar zi, dz, p, pd, h;
417 
418  zi = *(zgrj+i);
419  dz = z - zi;
420  if (std::abs(dz) < INTREPID_TOL) return 1.0;
421 
422  IntrepidPolylib::jacobfd (1, &zi, &p , (Scalar*)0, np-1, alpha+1, beta );
423  // need to use this routine in case z = -1 or 1
424  IntrepidPolylib::jacobd (1, &zi, &pd, np-1, alpha+1, beta );
425  h = (1.0 - zi)*pd - p;
426  IntrepidPolylib::jacobfd (1, &z, &p, (Scalar*)0, np-1, alpha+1, beta);
427  h = (1.0 - z )*p/(h*dz);
428 
429  return h;
430 }
431 
432 
433 template <class Scalar>
434 Scalar IntrepidPolylib::hglj (const int i, const Scalar z, const Scalar *zglj,
435  const int np, const Scalar alpha, const Scalar beta)
436 {
437  Scalar one = 1., two = 2.;
438  Scalar zi, dz, p, pd, h;
439 
440  zi = *(zglj+i);
441  dz = z - zi;
442  if (std::abs(dz) < INTREPID_TOL) return 1.0;
443 
444  IntrepidPolylib::jacobfd(1, &zi, &p , (Scalar*)0, np-2, alpha + one, beta + one);
445  // need to use this routine in case z = -1 or 1
446  IntrepidPolylib::jacobd (1, &zi, &pd, np-2, alpha + one, beta + one);
447  h = (one - zi*zi)*pd - two*zi*p;
448  IntrepidPolylib::jacobfd(1, &z, &p, (Scalar*)0, np-2, alpha + one, beta + one);
449  h = (one - z*z)*p/(h*dz);
450 
451  return h;
452 }
453 
454 
455 template <class Scalar>
456 void IntrepidPolylib::Imgj(Scalar *im, const Scalar *zgj, const Scalar *zm, const int nz,
457  const int mz, const Scalar alpha, const Scalar beta){
458  Scalar zp;
459  int i, j;
460 
461  /* old Polylib code */
462  for (i = 0; i < mz; ++i) {
463  zp = zm[i];
464  for (j = 0; j < nz; ++j) {
465  im[i*nz+j] = IntrepidPolylib::hgj(j, zp, zgj, nz, alpha, beta);
466  }
467  }
468 
469  /* original Nektar++ code: is this correct???
470  for (i = 0; i < nz; ++i) {
471  for (j = 0; j < mz; ++j)
472  {
473  zp = zm[j];
474  im [i*mz+j] = IntrepidPolylib::hgj(i, zp, zgj, nz, alpha, beta);
475  }
476  }
477  */
478 
479  return;
480 }
481 
482 
483 template <class Scalar>
484 void IntrepidPolylib::Imgrjm(Scalar *im, const Scalar *zgrj, const Scalar *zm, const int nz,
485  const int mz, const Scalar alpha, const Scalar beta)
486 {
487  Scalar zp;
488  int i, j;
489 
490  for (i = 0; i < mz; i++) {
491  zp = zm[i];
492  for (j = 0; j < nz; j++) {
493  im[i*nz+j] = IntrepidPolylib::hgrjm(j, zp, zgrj, nz, alpha, beta);
494  }
495  }
496 
497  /* original Nektar++ code: is this correct???
498  for (i = 0; i < nz; i++) {
499  for (j = 0; j < mz; j++)
500  {
501  zp = zm[j];
502  im [i*mz+j] = IntrepidPolylib::hgrjm(i, zp, zgrj, nz, alpha, beta);
503  }
504  }
505  */
506 
507  return;
508 }
509 
510 
511 template <class Scalar>
512 void IntrepidPolylib::Imgrjp(Scalar *im, const Scalar *zgrj, const Scalar *zm, const int nz,
513  const int mz, const Scalar alpha, const Scalar beta)
514 {
515  Scalar zp;
516  int i, j;
517 
518  for (i = 0; i < mz; i++) {
519  zp = zm[i];
520  for (j = 0; j < nz; j++) {
521  im [i*nz+j] = IntrepidPolylib::hgrjp(j, zp, zgrj, nz, alpha, beta);
522  }
523  }
524 
525  /* original Nektar++ code: is this correct?
526  for (i = 0; i < nz; i++) {
527  for (j = 0; j < mz; j++)
528  {
529  zp = zm[j];
530  im [i*mz+j] = IntrepidPolylib::hgrjp(i, zp, zgrj, nz, alpha, beta);
531  }
532  }
533  */
534 
535  return;
536 }
537 
538 
539 template <class Scalar>
540 void IntrepidPolylib::Imglj(Scalar *im, const Scalar *zglj, const Scalar *zm, const int nz,
541  const int mz, const Scalar alpha, const Scalar beta)
542 {
543  Scalar zp;
544  int i, j;
545 
546  for (i = 0; i < mz; i++) {
547  zp = zm[i];
548  for (j = 0; j < nz; j++) {
549  im[i*nz+j] = IntrepidPolylib::hglj(j, zp, zglj, nz, alpha, beta);
550  }
551  }
552 
553  /* original Nektar++ code: is this correct?
554  for (i = 0; i < nz; i++) {
555  for (j = 0; j < mz; j++)
556  {
557  zp = zm[j];
558  im[i*mz+j] = IntrepidPolylib::hglj(i, zp, zglj, nz, alpha, beta);
559  }
560  }
561  */
562 
563  return;
564 }
565 
566 
567 template <class Scalar>
568 void
570 jacobfd (const int np, const Scalar *z, Scalar *poly_in, Scalar *polyd,
571  const int n, const Scalar alpha, const Scalar beta)
572 {
573  const Scalar zero = 0.0, one = 1.0, two = 2.0;
574 
575  if (! np) {
576  return;
577  }
578 
579  if (n == 0) {
580  if (poly_in) {
581  for (int i = 0; i < np; ++i) {
582  poly_in[i] = one;
583  }
584  }
585  if (polyd) {
586  for (int i = 0; i < np; ++i) {
587  polyd[i] = zero;
588  }
589  }
590  }
591  else if (n == 1) {
592  if (poly_in) {
593  for (int i = 0; i < np; ++i) {
594  poly_in[i] = 0.5*(alpha - beta + (alpha + beta + two)*z[i]);
595  }
596  }
597  if (polyd) {
598  for (int i = 0; i < np; ++i) {
599  polyd[i] = 0.5*(alpha + beta + two);
600  }
601  }
602  }
603  else {
604  Scalar a1,a2,a3,a4;
605  Scalar apb = alpha + beta;
606  Scalar *poly, *polyn1,*polyn2;
607 
608  if (poly_in) { // switch for case of no poynomial function return
609  polyn1 = (Scalar *)malloc(2*np*sizeof(Scalar));
610  polyn2 = polyn1+np;
611  poly = poly_in;
612  }
613  else{
614  polyn1 = (Scalar *)malloc(3*np*sizeof(Scalar));
615  polyn2 = polyn1+np;
616  poly = polyn2+np;
617  }
618 
619  for (int i = 0; i < np; ++i) {
620  polyn2[i] = one;
621  polyn1[i] = 0.5*(alpha - beta + (alpha + beta + two)*z[i]);
622  }
623 
624  for (int k = 2; k <= n; ++k) {
625  a1 = two*k*(k + apb)*(two*k + apb - two);
626  a2 = (two*k + apb - one)*(alpha*alpha - beta*beta);
627  a3 = (two*k + apb - two)*(two*k + apb - one)*(two*k + apb);
628  a4 = two*(k + alpha - one)*(k + beta - one)*(two*k + apb);
629 
630  a2 /= a1;
631  a3 /= a1;
632  a4 /= a1;
633 
634  for (int i = 0; i < np; ++i) {
635  poly [i] = (a2 + a3*z[i])*polyn1[i] - a4*polyn2[i];
636  polyn2[i] = polyn1[i];
637  polyn1[i] = poly [i];
638  }
639  }
640 
641  if (polyd) {
642  a1 = n*(alpha - beta);
643  a2 = n*(two*n + alpha + beta);
644  a3 = two*(n + alpha)*(n + beta);
645  a4 = (two*n + alpha + beta);
646  a1 /= a4; a2 /= a4; a3 /= a4;
647 
648  // note polyn2 points to polyn1 at end of poly iterations
649  for (int i = 0; i < np; ++i) {
650  polyd[i] = (a1- a2*z[i])*poly[i] + a3*polyn2[i];
651  polyd[i] /= (one - z[i]*z[i]);
652  }
653  }
654 
655  free(polyn1);
656  }
657 }
658 
659 
660 template <class Scalar>
661 void IntrepidPolylib::jacobd(const int np, const Scalar *z, Scalar *polyd, const int n,
662  const Scalar alpha, const Scalar beta)
663 {
664  int i;
665  Scalar one = 1.0;
666  if(n == 0)
667  for(i = 0; i < np; ++i) polyd[i] = 0.0;
668  else{
669  //jacobf(np,z,polyd,n-1,alpha+one,beta+one);
670  IntrepidPolylib::jacobfd(np,z,polyd,(Scalar*)0,n-1,alpha+one,beta+one);
671  for(i = 0; i < np; ++i) polyd[i] *= 0.5*(alpha + beta + (Scalar)n + one);
672  }
673  return;
674 }
675 
676 
677 template <class Scalar>
678 void IntrepidPolylib::Jacobz(const int n, Scalar *z, const Scalar alpha, const Scalar beta){
679  int i,j,k;
680  Scalar dth = M_PI/(2.0*(Scalar)n);
681  Scalar poly,pder,rlast=0.0;
682  Scalar sum,delr,r;
683  Scalar one = 1.0, two = 2.0;
684 
685  if(!n)
686  return;
687 
688  for(k = 0; k < n; ++k){
689  r = -std::cos((two*(Scalar)k + one) * dth);
690  if(k) r = 0.5*(r + rlast);
691 
692  for(j = 1; j < INTREPID_POLYLIB_STOP; ++j){
693  IntrepidPolylib::jacobfd(1,&r,&poly, &pder, n, alpha, beta);
694 
695  for(i = 0, sum = 0.0; i < k; ++i) sum += one/(r - z[i]);
696 
697  delr = -poly / (pder - sum * poly);
698  r += delr;
699  if( std::abs(delr) < INTREPID_TOL ) break;
700  }
701  z[k] = r;
702  rlast = r;
703  }
704  return;
705 }
706 
707 
708 template <class Scalar>
709 void IntrepidPolylib::JacZeros(const int n, Scalar *a, const Scalar alpha, const Scalar beta){
710  int i;
711  Scalar apb, apbi,a2b2;
712  Scalar *b;
713 
714  if(!n)
715  return;
716 
717  b = (Scalar *) malloc(n*sizeof(Scalar));
718 
719  // generate normalised terms
720  apb = alpha + beta;
721  apbi = 2.0 + apb;
722 
723  b[n-1] = std::pow(2.0,apb+1.0)*IntrepidPolylib::gammaF(alpha+1.0)*IntrepidPolylib::gammaF(beta+1.0)/gammaF(apbi);
724  a[0] = (beta-alpha)/apbi;
725  b[0] = std::sqrt(4.0*(1.0+alpha)*(1.0+beta)/((apbi+1.0)*apbi*apbi));
726 
727  a2b2 = beta*beta-alpha*alpha;
728  for(i = 1; i < n-1; ++i){
729  apbi = 2.0*(i+1) + apb;
730  a[i] = a2b2/((apbi-2.0)*apbi);
731  b[i] = std::sqrt(4.0*(i+1)*(i+1+alpha)*(i+1+beta)*(i+1+apb)/
732  ((apbi*apbi-1)*apbi*apbi));
733  }
734 
735  apbi = 2.0*n + apb;
736  //a[n-1] = a2b2/((apbi-2.0)*apbi); // THIS IS A BUG!!!
737  if (n>1) a[n-1] = a2b2/((apbi-2.0)*apbi);
738 
739  // find eigenvalues
740  IntrepidPolylib::TriQL(n, a, b);
741 
742  free(b);
743  return;
744 }
745 
746 
747 template <class Scalar>
748 void IntrepidPolylib::TriQL(const int n, Scalar *d,Scalar *e) {
749  int m,l,iter,i,k;
750  Scalar s,r,p,g,f,dd,c,b;
751 
752  for (l=0;l<n;l++) {
753  iter=0;
754  do {
755  for (m=l;m<n-1;m++) {
756  dd=std::abs(d[m])+std::abs(d[m+1]);
757  if (std::abs(e[m])+dd == dd) break;
758  }
759  if (m != l) {
760  if (iter++ == 50){
761  TEUCHOS_TEST_FOR_EXCEPTION((1),
762  std::runtime_error,
763  ">>> ERROR (IntrepidPolylib): Too many iterations in TQLI.");
764  }
765  g=(d[l+1]-d[l])/(2.0*e[l]);
766  r=std::sqrt((g*g)+1.0);
767  //g=d[m]-d[l]+e[l]/(g+sign(r,g));
768  g=d[m]-d[l]+e[l]/(g+((g)<0 ? -std::abs(r) : std::abs(r)));
769  s=c=1.0;
770  p=0.0;
771  for (i=m-1;i>=l;i--) {
772  f=s*e[i];
773  b=c*e[i];
774  if (std::abs(f) >= std::abs(g)) {
775  c=g/f;
776  r=std::sqrt((c*c)+1.0);
777  e[i+1]=f*r;
778  c *= (s=1.0/r);
779  } else {
780  s=f/g;
781  r=std::sqrt((s*s)+1.0);
782  e[i+1]=g*r;
783  s *= (c=1.0/r);
784  }
785  g=d[i+1]-p;
786  r=(d[i]-g)*s+2.0*c*b;
787  p=s*r;
788  d[i+1]=g+p;
789  g=c*r-b;
790  }
791  d[l]=d[l]-p;
792  e[l]=g;
793  e[m]=0.0;
794  }
795  } while (m != l);
796  }
797 
798  // order eigenvalues
799  for(i = 0; i < n-1; ++i){
800  k = i;
801  p = d[i];
802  for(l = i+1; l < n; ++l)
803  if (d[l] < p) {
804  k = l;
805  p = d[l];
806  }
807  d[k] = d[i];
808  d[i] = p;
809  }
810 }
811 
812 
813 template <class Scalar>
814 Scalar IntrepidPolylib::gammaF(const Scalar x){
815  Scalar gamma = 1.0;
816 
817  if (x == -0.5) gamma = -2.0*std::sqrt(M_PI);
818  else if (!x) return gamma;
819  else if ((x-(int)x) == 0.5){
820  int n = (int) x;
821  Scalar tmp = x;
822 
823  gamma = std::sqrt(M_PI);
824  while(n--){
825  tmp -= 1.0;
826  gamma *= tmp;
827  }
828  }
829  else if ((x-(int)x) == 0.0){
830  int n = (int) x;
831  Scalar tmp = x;
832 
833  while(--n){
834  tmp -= 1.0;
835  gamma *= tmp;
836  }
837  }
838  else
839  TEUCHOS_TEST_FOR_EXCEPTION((1),
840  std::invalid_argument,
841  ">>> ERROR (IntrepidPolylib): Argument is not of integer or half order.");
842  return gamma;
843 }
844 
845 } // end of namespace Intrepid
846 
static void Dgrjm(Scalar *D, const Scalar *z, const int np, const Scalar alpha, const Scalar beta)
Compute the Derivative Matrix and its transpose associated with the Gauss-Radau-Jacobi zeros with a z...
static void Imgj(Scalar *im, const Scalar *zgj, const Scalar *zm, const int nz, const int mz, const Scalar alpha, const Scalar beta)
Interpolation Operator from Gauss-Jacobi points to an arbitrary distribution at points zm...
static Scalar hgrjp(const int i, const Scalar z, const Scalar *zgrj, const int np, const Scalar alpha, const Scalar beta)
Compute the value of the i th Lagrangian interpolant through the np Gauss-Radau-Jacobi points zgrj at...
static void zwgrjp(Scalar *z, Scalar *w, const int np, const Scalar alpha, const Scalar beta)
Gauss-Radau-Jacobi zeros and weights with end point at z=1.
static void TriQL(const int n, Scalar *d, Scalar *e)
QL algorithm for symmetric tridiagonal matrix.
#define INTREPID_POLYLIB_STOP
Maximum number of iterations in polynomial defalation routine Jacobz.
static Scalar hglj(const int i, const Scalar z, const Scalar *zglj, const int np, const Scalar alpha, const Scalar beta)
Compute the value of the i th Lagrangian interpolant through the np Gauss-Lobatto-Jacobi points zglj ...
static void zwglj(Scalar *z, Scalar *w, const int np, const Scalar alpha, const Scalar beta)
Gauss-Lobatto-Jacobi zeros and weights with end point at z=-1,1.
static Scalar hgj(const int i, const Scalar z, const Scalar *zgj, const int np, const Scalar alpha, const Scalar beta)
Compute the value of the i th Lagrangian interpolant through the np Gauss-Jacobi points zgj at the ar...
static void zwgrjm(Scalar *z, Scalar *w, const int np, const Scalar alpha, const Scalar beta)
Gauss-Radau-Jacobi zeros and weights with end point at z=-1.
static Scalar hgrjm(const int i, const Scalar z, const Scalar *zgrj, const int np, const Scalar alpha, const Scalar beta)
Compute the value of the i th Lagrangian interpolant through the np Gauss-Radau-Jacobi points zgrj at...
static void JacZeros(const int n, Scalar *a, const Scalar alpha, const Scalar beta)
Zero determination through the eigenvalues of a tridiagonal matrix from the three term recursion rela...
static void jacobd(const int np, const Scalar *z, Scalar *polyd, const int n, const Scalar alpha, const Scalar beta)
Calculate the derivative of Jacobi polynomials.
static void Dgrjp(Scalar *D, const Scalar *z, const int np, const Scalar alpha, const Scalar beta)
Compute the Derivative Matrix associated with the Gauss-Radau-Jacobi zeros with a zero at z=1...
static void Imgrjm(Scalar *im, const Scalar *zgrj, const Scalar *zm, const int nz, const int mz, const Scalar alpha, const Scalar beta)
Interpolation Operator from Gauss-Radau-Jacobi points (including z=-1) to an arbitrary distrubtion at...
static void jacobfd(const int np, const Scalar *z, Scalar *poly_in, Scalar *polyd, const int n, const Scalar alpha, const Scalar beta)
Routine to calculate Jacobi polynomials, , and their first derivative, .
static void Jacobz(const int n, Scalar *z, const Scalar alpha, const Scalar beta)
Calculate the n zeros, z, of the Jacobi polynomial, i.e. .
static void Dgj(Scalar *D, const Scalar *z, const int np, const Scalar alpha, const Scalar beta)
Compute the Derivative Matrix and its transpose associated with the Gauss-Jacobi zeros.
static void Dglj(Scalar *D, const Scalar *z, const int np, const Scalar alpha, const Scalar beta)
Compute the Derivative Matrix associated with the Gauss-Lobatto-Jacobi zeros.
static void zwgj(Scalar *z, Scalar *w, const int np, const Scalar alpha, const Scalar beta)
Gauss-Jacobi zeros and weights.
static void Imglj(Scalar *im, const Scalar *zglj, const Scalar *zm, const int nz, const int mz, const Scalar alpha, const Scalar beta)
Interpolation Operator from Gauss-Lobatto-Jacobi points to an arbitrary distrubtion at points zm...
static Scalar gammaF(const Scalar x)
Calculate the Gamma function , , for integer values x and halves.
#define jacobz(n, z, alpha, beta)
zero determination using Newton iteration with polynomial deflation
static void Imgrjp(Scalar *im, const Scalar *zgrj, const Scalar *zm, const int nz, const int mz, const Scalar alpha, const Scalar beta)
Interpolation Operator from Gauss-Radau-Jacobi points (including z=1) to an arbitrary distrubtion at ...