SundanceFeketeTriangleQuadrature.cpp
Go to the documentation of this file.
00001 /* @HEADER@ */
00002 // ************************************************************************
00003 // 
00004 //                             Sundance
00005 //                 Copyright 2011 Sandia Corporation
00006 // 
00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00008 // the U.S. Government retains certain rights in this software.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Kevin Long (kevin.long@ttu.edu)
00038 // 
00039 
00040 /* @HEADER@ */
00041 
00042 
00043 #include "SundanceFeketeTriangleQuadrature.hpp"
00044 #include "PlayaExceptions.hpp"
00045 #include "SundanceOut.hpp"
00046 #include "SundancePoint.hpp"
00047 #include "PlayaTabs.hpp"
00048 
00049 using namespace Sundance;
00050 using namespace Teuchos;
00051 
00052 extern "C"
00053 {
00054 /* matrix-vector multiplication */
00055 void dgemv_(char *trans, int *m, int *n, double *alpha, double *a, int *lda,
00056     double *x, int *incx, double *beta, double *y, int *incy);
00057 
00058 /* LAPACK factorization */
00059 void dgetrf_(const int* M, const int* N, double* A, const int* lda,
00060     const int* iPiv, int* info);
00061 
00062 /* LAPACK inversion of factorized matrix */
00063 void dgetri_(const int* n, double* a, const int* lda, const int* iPiv,
00064     double* work, const int* lwork, int* info);
00065 }
00066 
00067 /**
00068  *  Reference:
00069  *  T. Warburton, An explicit construction of interpolation nodes on the simplex
00070  *  J. Eng. Math. (2006) 56, pp. 247-262
00071  */
00072 void FeketeTriangleQuadrature::getPoints(int order, Array<double>& wgt, Array<
00073     double>& x, Array<double>& y)
00074 {
00075   Array<double> w;
00076   Array<int> multiplicity;
00077   Array<Array<double> > q;
00078 
00079   if (order == 1)
00080   {
00081     // Delete it? One gets 2nd order for the same price...
00082     multiplicity = tuple(3);
00083     q.resize(1);
00084     w = tuple(1.0 / 3.0);
00085     q[0] = tuple(1.0, 0.0, 0.0);
00086   }
00087   else if (order == 2)
00088   {
00089     // Corners have weights == 0, but we list them anyway for adaptive cell integration
00090     multiplicity = tuple(3, 3);
00091     q.resize(2);
00092     w = tuple(0.0, 1.0 / 3.0);
00093     q[0] = tuple(1.0, 0.0, 0.0);
00094     q[1] = tuple(0.0, 0.5, 0.5);
00095   }
00096   else if (order == 3)
00097   {
00098     multiplicity = tuple(1, 3, 6);
00099     q.resize(3);
00100     w = tuple(0.45, 1.0 / 60.0, 1.0 / 12.0);
00101     q[0] = tuple(1.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0);
00102     q[1] = tuple(1.0, 0.0, 0.0);
00103     q[2] = tuple(0.000000000000000, 0.276393202250021, 0.723606797749979);
00104   }
00105   else if (order == 4)
00106   {
00107     multiplicity = tuple(3, 3, 3, 6);
00108     q.resize(4);
00109     w = tuple(-0.002443433685378, 0.040684938686721, 0.200360939516545,
00110         0.047365444407723);
00111     q[0] = tuple(1.0, 0.0, 0.0);
00112     q[1] = tuple(0.0, 0.5, 0.5);
00113     q[2] = tuple(0.551583507555305, 0.224208246222347, 0.224208246222347);
00114     q[3] = tuple(0.000000000000000, 0.172673164646012, 0.827326835353988);
00115   }
00116   else if (order == 5)
00117   {
00118     multiplicity = tuple(3, 3, 3, 6, 6);
00119     q.resize(5);
00120     w = tuple(0.005695992854379, 0.125094845896272, 0.116460019007601,
00121         0.012270695896559, 0.030770541890982);
00122     q[0] = tuple(1.0, 0.0, 0.0);
00123     q[1] = tuple(0.684472514501908, 0.157763742749046, 0.157763742749046);
00124     q[2] = tuple(0.171245477332074, 0.414377261333963, 0.414377261333963);
00125     q[3] = tuple(0.000000000000000, 0.117472338035268, 0.882527661964732);
00126     q[4] = tuple(0.000000000000000, 0.357384241759678, 0.642615758240322);
00127   }
00128   else if (order == 6)
00129   {
00130     multiplicity = tuple(1, 3, 3, 3, 6, 6, 6);
00131     q.resize(7);
00132     w = tuple(0.078877036160935, -0.001814501354491, 0.021979435329947,
00133         0.055816419204784, 0.013312944349338, 0.013197985920130,
00134         0.089018887113590);
00135     q[0] = tuple(1.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0);
00136     q[1] = tuple(1.0, 0.0, 0.0);
00137     q[2] = tuple(0.0, 0.5, 0.5);
00138     q[3] = tuple(0.769520861253251, 0.115239569373374, 0.115239569373375);
00139     q[4] = tuple(0.000000000000000, 0.084888051860717, 0.915111948139283);
00140     q[5] = tuple(0.000000000000000, 0.265575603264643, 0.734424396735357);
00141     q[6] = tuple(0.127944523240232, 0.317617605129902, 0.554437871629866);
00142   }
00143   else if (order == 9)
00144   {
00145     multiplicity = tuple(1, 3, 3, 3, 3, 6, 6, 6, 6, 6, 6, 6);
00146     q.resize(12);
00147     w = tuple(0.054830037550851, 0.001827879484114, 0.019001036234018,
00148         0.036420799722133, 0.030554508665561, -0.000151089387317,
00149         0.005590374483628, 0.004964725734335, 0.006760863782205,
00150         0.020422905832759, 0.035099032364350, 0.040939402211986);
00151     q[0] = tuple(1.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0);
00152     q[1] = tuple(1.0, 0.0, 0.0);
00153     q[2] = tuple(0.888511356254118, 0.055744321872941, 0.055744321872941);
00154     q[3] = tuple(0.643273196352075, 0.178363401823962, 0.178363401823962);
00155     q[4] = tuple(0.067157296821421, 0.466421351589290, 0.466421351589290);
00156     q[5] = tuple(0.000000000000000, 0.040233045916771, 0.959766954083229);
00157     q[6] = tuple(0.000000000000000, 0.130613067447248, 0.869386932552752);
00158     q[7] = tuple(0.000000000000000, 0.261037525094778, 0.738962474905222);
00159     q[8] = tuple(0.000000000000000, 0.417360521166807, 0.582639478833193);
00160     q[9] = tuple(0.063288094133481, 0.162084689378845, 0.774627216487674);
00161     q[10] = tuple(0.066372381175690, 0.304692583561054, 0.628935035263256);
00162     q[11] = tuple(0.185067879630320, 0.326702931348863, 0.488229189020817);
00163   }
00164   /* Points and weights according to Fekete approach by
00165    * Taylor, Wingate, Vincent 2000
00166    * Unfortunately not sufficient digits for type 'double'!
00167    else if (order == 6)
00168    {
00169    multiplicity = tuple(1, 3, 3, 3, 6, 6, 6);
00170    q.resize(7);
00171    w = tuple(0.2178563571, 0.1104193374, 0.0358939762, 0.0004021278,
00172    0.1771348660, 0.0272344079, 0.0192969460);
00173    q[0] = tuple(0.3333333333, 0.3333333333, 0.3333333334);
00174    q[1] = tuple(0.7873290632, 0.1063354684, 0.1063354684);
00175    q[2] = tuple(0.0000000000, 0.5000000000, 0.5000000000);
00176    q[3] = tuple(1.0000000000, 0.0000000000, 0.0000000000);
00177    q[4] = tuple(0.1171809171, 0.3162697959, 0.5665492870);
00178    q[5] = tuple(0.0000000000, 0.2655651402, 0.7344348598);
00179    q[6] = tuple(0.0000000000, 0.0848854223, 0.9151145777);
00180    }
00181    else if (order == 9)
00182    {
00183    multiplicity = tuple(1, 3, 3, 3, 3, 6, 6, 6, 6, 6, 6, 6);
00184    q.resize(12);
00185    w = tuple(0.1096011288, 0.0767491008, 0.0646677819, 0.0276211659,
00186    0.0013925011, 0.0933486453, 0.0619010169, 0.0437466450,
00187    0.0114553907, 0.0093115568, 0.0078421987, 0.0022457501);
00188    q[0] = tuple(0.3333333333, 0.3333333333, 0.3333333334);
00189    q[1] = tuple(0.6591363598, 0.1704318201, 0.1704318201);
00190    q[2] = tuple(0.0600824712, 0.4699587644, 0.4699587644);
00191    q[3] = tuple(0.9021308608, 0.0489345696, 0.0489345696);
00192    q[4] = tuple(1.0000000000, 0.0000000000, 0.0000000000);
00193    q[5] = tuple(0.1784337588, 0.3252434900, 0.4963227512);
00194    q[6] = tuple(0.0588564879, 0.3010242110, 0.6401193011);
00195    q[7] = tuple(0.0551758079, 0.1543901944, 0.7904339977);
00196    q[8] = tuple(0.0000000000, 0.4173602935, 0.5826397065);
00197    q[9] = tuple(0.0000000000, 0.2610371960, 0.7389628040);
00198    q[10] = tuple(0.0000000000, 0.1306129092, 0.8693870908);
00199    q[11] = tuple(0.0000000000, 0.0402330070, 0.9597669930);
00200    }*/
00201 
00202   else
00203   {
00204 #ifndef TRILINOS_7
00205     SUNDANCE_ERROR("symmetric Fekete quadrature rule order "
00206         << order <<
00207         " for triangles not available");
00208 #else
00209     SUNDANCE_ERROR7("symmetric Fekete quadrature rule order "
00210         << order <<
00211         " for triangles not available");
00212 #endif
00213   }
00214 
00215   for (int i = 0; i < q.length(); i++)
00216   {
00217     Array<Array<double> > qPerm;
00218     permute(multiplicity[i], q[i], qPerm);
00219     for (int j = 0; j < multiplicity[i]; j++)
00220     {
00221       x.append(qPerm[j][0]);
00222       y.append(qPerm[j][1]);
00223       wgt.append(w[i]);
00224     }
00225   }
00226 
00227 }
00228 
00229 bool FeketeTriangleQuadrature::supportsOrder(int order)
00230 {
00231   if ((order >= 1 && order <= 6) || order == 9)
00232     return true;
00233   return false;
00234 }
00235 
00236 /**
00237  * Here we calculate coefficients for Proriol-Koornwinder-Dubiner polynomials
00238  * so that they form a Lagrange basis at given (Fekete quadrature) points in the
00239  * triangle
00240  */
00241 void FeketeTriangleQuadrature::computeBasisCoeffs(const int order, Array<double>& basisCoeffs)
00242 {
00243   // Get Fekete points of chosen order
00244   Array<Point> feketePts;
00245 
00246   Array<double> x;
00247   Array<double> y;
00248   Array<double> w;
00249   getPoints(order, w, x, y);
00250   int nFeketePts = w.length();
00251   feketePts.resize(nFeketePts);
00252   for (int i = 0; i < nFeketePts; i++)
00253     feketePts[i] = Point(x[i], y[i]);
00254 
00255   // We construct a Lagrange basis at feketePts (there are nFeketePts of them)
00256   // Each basis polynomial itself is given by a linear combination of
00257   // nFeketePts PKD polynomials and their coefficients
00258   basisCoeffs.resize(nFeketePts * nFeketePts);
00259 
00260   // Let's compute the coefficients:
00261   // Build Vandermonde matrix of PKD basis at Fekete points
00262   for (int n = 0; n < nFeketePts; n++)
00263   {
00264     // Set pointer to beginning of n-th row and determine values of
00265     // PKD polynomials at n-th Fekete point
00266     double* start = &(basisCoeffs[n * nFeketePts]);
00267     evalPKDpolynomials(order, feketePts[n][0], feketePts[n][1], start);
00268   }
00269 
00270   // Invert Vandermonde matrix to obtain basis coefficients:
00271   // LAPACK error flag, array for switched rows, work array
00272   int lapack_err = 0;
00273   Array<int> pivot;
00274   pivot.resize(nFeketePts);
00275   Array<double> work;
00276   work.resize(1);
00277   int lwork = -1;
00278 
00279   // LU factorization
00280   ::dgetrf_(&nFeketePts, &nFeketePts, &(basisCoeffs[0]), &nFeketePts,
00281       &(pivot[0]), &lapack_err);
00282 
00283   TEUCHOS_TEST_FOR_EXCEPTION(
00284       lapack_err != 0,
00285       std::runtime_error,
00286       "FeketeTriangleQuadrature::computeBasisCoeffs(): factorization of generalized Vandermonde matrix failed");
00287 
00288   // Determine work array size and invert factorized matrix
00289   ::dgetri_(&nFeketePts, &(basisCoeffs[0]), &nFeketePts, &(pivot[0]),
00290       &(work[0]), &lwork, &lapack_err);
00291   lwork = (int) work[0];
00292   work.resize(lwork);
00293   ::dgetri_(&nFeketePts, &(basisCoeffs[0]), &nFeketePts, &(pivot[0]),
00294       &(work[0]), &lwork, &lapack_err);
00295 
00296   TEUCHOS_TEST_FOR_EXCEPTION(
00297       lapack_err != 0,
00298       std::runtime_error,
00299       "FeketeTriangleQuadrature::computeBasisCoeffs(): inversion of generalized Vandermonde matrix failed");
00300 }
00301 
00302 /**
00303  * Evaluates all basis functions of a Proriol-Koornwinder-Dubiner basis
00304  * up to the given order at (x,y) in reference (barycentric) coordinates of a triangle;
00305  * Missing third coordinate z = 1-x-y
00306  */
00307 void FeketeTriangleQuadrature::evalPKDpolynomials(int order, double x,
00308     double y, double* resultPtr)
00309 {
00310   int i = 0;
00311   double Legendre = 1.0;
00312   double Legendre1 = 0.0;
00313   for (int k = 0; k <= order; k++)
00314   {
00315     double Jacobi = 1.0;
00316     double Jacobi1 = 0.0;
00317 
00318     // 'i' implements a bijection polynomial index (k,l) onto basis index (i)
00319     for (int l = 0; l <= order - k; l++)
00320     {
00321       // Overall we will write (order+1)*(order+2)/2 values
00322       resultPtr[i++] = Legendre * pow(x + y, double(k)) * Jacobi;
00323 
00324       // Update Jacobi polynomial
00325       double Jacobi2 = Jacobi1;
00326       Jacobi1 = Jacobi;
00327       int c1 = 2 * (l + 1) * (l + 2 * k + 2) * (2 * l + 2 * k + 1);
00328       int c2 = (2 * l + 2 * k + 2) * (2 * k + 1) * (2 * k + 1);
00329       int c3 = (2 * l + 2 * k + 1) * (2 * l + 2 * k + 2) * (2 * l + 2 * k
00330           + 3);
00331       int c4 = 2 * (l + 2 * k + 1) * l * (2 * l + 2 * k + 3);
00332       Jacobi = ((c2 + c3 * (1.0 - 2.0 * x - 2.0 * y)) * Jacobi1 - c4
00333           * Jacobi2) / c1;
00334     }
00335 
00336     // Update Legendre polynomial
00337     double Legendre2 = Legendre1;
00338     Legendre1 = Legendre;
00339     double xy = x + y;
00340     if (::fabs(xy) > 0.0)
00341     {
00342       xy = (2.0 * k + 1.0) * ((x - y) / xy) * Legendre1;
00343     }
00344     Legendre = (xy - k * Legendre2) / (k + 1);
00345   }
00346 }
00347 
00348 void FeketeTriangleQuadrature::permute(int m, const Array<double>& q, Array<
00349     Array<double> >& qPerm)
00350 {
00351   qPerm.resize(m);
00352   if (m == 1)
00353   {
00354     qPerm[0] = q;
00355   }
00356   else if (m == 3)
00357   {
00358     qPerm[0] = tuple(q[0], q[1], q[2]);
00359     qPerm[1] = tuple(q[1], q[0], q[2]);
00360     qPerm[2] = tuple(q[2], q[1], q[0]);
00361   }
00362   else if (m == 6)
00363   {
00364     qPerm[0] = tuple(q[0], q[1], q[2]);
00365     qPerm[1] = tuple(q[0], q[2], q[1]);
00366     qPerm[2] = tuple(q[1], q[0], q[2]);
00367     qPerm[3] = tuple(q[1], q[2], q[0]);
00368     qPerm[4] = tuple(q[2], q[1], q[0]);
00369     qPerm[5] = tuple(q[2], q[0], q[1]);
00370   }
00371   else
00372   {
00373 #ifndef TRILINOS_7
00374     SUNDANCE_ERROR("invalid multiplicity "
00375         << m <<
00376         " in FeketeTriangleQuadrature::permute()");
00377 #else
00378     SUNDANCE_ERROR7("invalid multiplicity "
00379         << m <<
00380         " in FeketeTriangleQuadrature::permute()");
00381 #endif
00382   }
00383 }
00384 
00385 bool FeketeTriangleQuadrature::test(int p)
00386 {
00387   Array<double> w;
00388   Array<double> x;
00389   Array<double> y;
00390 
00391   getPoints(p, w, x, y);
00392   bool pass = true;
00393 
00394   for (int a = 0; a <= p; a++)
00395   {
00396     for (int b = 0; b < p - a; b++)
00397     {
00398       int cMax = p - a - b;
00399       for (int c = 0; c <= cMax; c++)
00400       {
00401         double sum = 0.0;
00402         for (int q = 0; q < w.length(); q++)
00403         {
00404           sum += 0.5 * w[q] * pow(x[q], (double) a) * pow(y[q],
00405               (double) b) * pow(1.0 - x[q] - y[q], (double) c);
00406         }
00407         double err = fabs(sum - exact(a, b, c));
00408         bool localPass = err < 1.0e-14;
00409         pass = pass && localPass;
00410         if (!localPass)
00411         {
00412           fprintf(
00413               stderr,
00414               "order=%d m (%d, %d, %d) q=%22.15g exact=%22.15g\n",
00415               p, a, b, c, sum, exact(a, b, c));
00416           std::cerr << "error = " << err << std::endl;
00417         }
00418       }
00419     }
00420   }
00421   return pass;
00422 }
00423 
00424 double FeketeTriangleQuadrature::exact(int a, int b, int c)
00425 {
00426   return fact(a) * fact(b) * fact(c) / fact(a + b + c + 2);
00427 }
00428 
00429 double FeketeTriangleQuadrature::fact(int x)
00430 {
00431   if (x == 0)
00432     return 1.0;
00433   return x * fact(x - 1);
00434 }
00435 

Site Contact