SundanceTriangleQuadrature.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 "SundanceTriangleQuadrature.hpp"
00044 #include "PlayaExceptions.hpp"
00045 #include "SundanceOut.hpp"
00046 #include "SundanceGauss1D.hpp"
00047 #include "PlayaTabs.hpp"
00048 
00049 using namespace Sundance;
00050 using namespace Sundance;
00051 using namespace Sundance;
00052 using namespace Sundance;
00053 using namespace Sundance;
00054 using namespace Teuchos;
00055 
00056 void TriangleQuadrature::getPoints(int order, Array<double>& wgt,
00057                                    Array<double>& x,
00058                                    Array<double>& y)
00059 {
00060   if (!getSymmetricPoints(order, wgt, x, y)) 
00061     {
00062       getNonsymmetricPoints(order, wgt, x, y);
00063     }
00064 }
00065 
00066 
00067 bool TriangleQuadrature::getSymmetricPoints(int order, Array<double>& wgt,
00068                                             Array<double>& x,
00069                                             Array<double>& y)
00070 {
00071 
00072   int np;
00073   Array<double> w;
00074   Array<int> multiplicity;
00075   Array<Array<double> > q;
00076 
00077   if (order==1)
00078     {
00079       multiplicity = tuple(1);
00080       np = 1;
00081       w = tuple(1.0);
00082       q.resize(1);
00083       q[0] = tuple(1.0/3.0, 1.0/3.0, 1.0/3.0);
00084     }
00085   else if (order==2)
00086     {
00087       multiplicity = tuple(3);
00088       np = 3;
00089       w = tuple(1.0/3.0);
00090       q.resize(1);
00091       q[0] = tuple(2.0/3.0, 1.0/6.0, 1.0/6.0);
00092     }
00093   else if (order==3)
00094     {
00095       multiplicity = tuple(6);
00096       np = 6;
00097       w = tuple(1.0/6.0);
00098       q.resize(1);
00099       q[0] = tuple(0.659027622374092, 0.231933368553031, 0.109039009072877);
00100     }
00101   else if (order==4)
00102     {
00103       multiplicity = tuple(3, 3);
00104       np = 6;
00105       w = tuple(0.109951743655322, 0.223381589678011);
00106       q.resize(2);
00107       q[0] = tuple(0.816847572980459, 0.091576213509771, 0.091576213509771);
00108       q[1] = tuple(0.108103018168070, 0.445948490915965, 0.445948490915965);
00109     }
00110   else if (order==5)
00111     {
00112       multiplicity = tuple(1, 3, 3);
00113       np = 7;
00114       q.resize(3);
00115       w = tuple(0.22500000000000, 0.125939180544827, 0.132394152788506);
00116       q[0] = tuple(1.0/3.0, 1.0/3.0, 1.0/3.0);
00117       q[1] = tuple(0.797426985353087, 0.101286507323456, 0.101286507323456);
00118       q[2] = tuple(0.059715871789770, 0.470142064105115, 0.470142064105115);
00119     }
00120   else if (order==6)
00121     {
00122       multiplicity = tuple(3, 3, 6);
00123       np = 12;
00124       q.resize(3);
00125       w = tuple(0.050844906370207, 0.116786275726379, 0.082851075618374);
00126       q[0] = tuple(0.873821971016996, 0.063089014491502, 0.063089014491502);
00127       q[1] = tuple(0.501426509658179, 0.249286745170910, 0.249286745170910);
00128       q[2] = tuple(0.636502499121399, 0.310352451033784, 0.053145049844817);
00129     }
00130   else
00131     {
00132       return false;
00133     }
00134 
00135   for (int i=0; i<q.length(); i++)
00136     {
00137       Array<Array<double> > qPerm;
00138       permute(multiplicity[i], q[i], qPerm);
00139       for (int j=0; j<multiplicity[i]; j++)
00140         {
00141           x.append(qPerm[j][0]);
00142           y.append(qPerm[j][1]);
00143           wgt.append(w[i]);
00144         }
00145     }
00146   
00147   return true;
00148 }
00149 
00150 
00151 
00152 
00153 void TriangleQuadrature::getNonsymmetricPoints(int order, Array<double>& wgt,
00154                                                Array<double>& x,
00155                                                Array<double>& y)
00156 {
00157   int nNodes = (order+3)/2;
00158   Gauss1D rule(nNodes, -1.0, 1.0);
00159   Array<double> s = rule.nodes();
00160   Array<double> t = s;
00161   Array<double> w = rule.weights();
00162   int n = rule.nPoints();
00163 
00164   wgt.resize(n*n);
00165   x.resize(n*n);
00166   y.resize(n*n);
00167 
00168   int k=0;
00169   for (int i=0; i<n; i++)
00170     {
00171       double p = (1.0+s[i])/2.0;
00172       double J = 1.0-p;
00173       for (int j=0; j<n; j++, k++)
00174         {
00175           double q = (1.0 - p)*(1.0+t[j])/2.0;
00176           x[k] = p;
00177           y[k] = q;
00178           wgt[k] = 0.5*w[i]*w[j]*J;
00179         }
00180     }
00181 }
00182 
00183 
00184 void TriangleQuadrature::permute(int m, const Array<double>& q,
00185                                  Array<Array<double> >& qPerm)
00186 {
00187   qPerm.resize(m);
00188   if (m==1)
00189     {
00190       qPerm[0] = q;
00191     }
00192   else if (m==3)
00193     {
00194       qPerm[0] = tuple(q[0], q[1], q[2]);
00195       qPerm[1] = tuple(q[1], q[0], q[2]);
00196       qPerm[2] = tuple(q[2], q[1], q[0]);
00197     }
00198   else if (m==6)
00199     {
00200       qPerm[0] = tuple(q[0], q[1], q[2]);
00201       qPerm[1] = tuple(q[0], q[2], q[1]);
00202       qPerm[2] = tuple(q[1], q[0], q[2]);
00203       qPerm[3] = tuple(q[1], q[2], q[0]);
00204       qPerm[4] = tuple(q[2], q[1], q[0]);
00205       qPerm[5] = tuple(q[2], q[0], q[1]);
00206     }
00207   else
00208     {
00209 #ifndef TRILINOS_7
00210       SUNDANCE_ERROR("invalid multiplicity " 
00211                      << m <<
00212                      " in TriangleQuadrature::permute()");
00213 #else
00214       SUNDANCE_ERROR7("invalid multiplicity " 
00215                      << m <<
00216                      " in TriangleQuadrature::permute()");
00217 #endif
00218     }
00219 }
00220 
00221 bool TriangleQuadrature::test(int p)
00222 {
00223   Array<double> w;
00224   Array<double> x;
00225   Array<double> y;
00226 
00227   getPoints(p, w, x, y);
00228   bool pass = true;
00229   
00230   for (int a=0; a<=p; a++)
00231     {
00232       for (int b=0; b<p-a; b++)
00233         {
00234           int cMax = p - a - b;
00235           for (int c=0; c<=cMax; c++)
00236             {
00237               double sum = 0.0;
00238               for (int q=0; q<w.length(); q++)
00239                 {
00240                   sum += 0.5*w[q] * pow(x[q], (double) a) * pow(y[q], (double) b)
00241                     * pow(1.0 - x[q] - y[q], (double) c);
00242                 }
00243               double err = fabs(sum - exact(a,b,c));
00244               bool localPass = err < 1.0e-14;
00245               pass = pass && localPass;
00246               if (!localPass)
00247                 {
00248                   fprintf(stderr, "order=%d m (%d, %d, %d) q=%22.15g exact=%22.15g\n", p, a, b, c, sum, exact(a, b, c));
00249                   std::cerr << "error = " << err << std::endl;
00250                 }
00251             }
00252         }
00253     }
00254   return pass;
00255 }
00256 
00257 double TriangleQuadrature::exact(int a, int b, int c)
00258 {
00259   return fact(a)*fact(b)*fact(c)/fact(a+b+c+2);
00260 }
00261 
00262 double TriangleQuadrature::fact(int x)
00263 {
00264   if (x==0) return 1.0;
00265   return x*fact(x-1);
00266 }
00267 
00268 

Site Contact