SundanceTetQuadrature.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 "SundanceTetQuadrature.hpp"
00044 #include "PlayaExceptions.hpp"
00045 #include "SundanceOut.hpp"
00046 
00047 using namespace Sundance;
00048 using namespace Sundance;
00049 using namespace Sundance;
00050 using namespace Sundance;
00051 using namespace Sundance;
00052 using namespace Teuchos;
00053 
00054 
00055 void TetQuadrature::getPoints(int order, Array<double>& wgt,
00056                               Array<double>& x,
00057                               Array<double>& y,
00058                               Array<double>& z)
00059 {
00060   Array<double> w;
00061   Array<int> multiplicity;
00062   Array<Array<double> > q;
00063 
00064   if (order==1)
00065     {
00066       multiplicity = tuple(1);
00067       w = tuple(1.0);
00068       q.resize(1);
00069       q[0] = tuple(0.25);
00070     }
00071   else if (order==2)
00072     {
00073       multiplicity = tuple(4);
00074       w = tuple(0.25);
00075       q.resize(1);
00076       q[0] = tuple(0.5854101966249685, 0.1381966011250105);
00077     }
00078   else if (order==4)
00079     {
00080       multiplicity = tuple(4, 12);
00081       w = tuple(0.05037379410012282, 0.06654206863329239);
00082       q.resize(2);
00083       q[0] = tuple(0.7716429020672371, 0.7611903264425430e-01);
00084       q[1] = tuple(0.1197005277978019, 0.7183164526766925e-01, 0.4042339134672644);
00085     }
00086   else if (order==6)
00087     {
00088       multiplicity = tuple(1, 4, 12, 12);
00089       w = tuple(0.9040129046014750e-01, 0.1911983427899124e-01,
00090                 0.4361493840666568e-01, 0.2581167596199161e-01);
00091       q.resize(4);
00092       q[0] = tuple(0.25);
00093       q[1] = tuple(0.8277192480479295, 0.5742691731735683e-01);
00094       q[2] = tuple(0.5135188412556341e-01, 0.4860510285706072, 0.2312985436519147);
00095       q[3] = tuple(0.2967538129690260, 0.6081079894015281, 0.4756909881472290e-01);
00096     }
00097   else
00098     {
00099 #ifndef TRILINOS_7
00100       SUNDANCE_ERROR("symmetric quadrature rule order " 
00101                      << order << 
00102                      " not available for triangles");
00103 #else
00104       SUNDANCE_ERROR7("symmetric quadrature rule order " 
00105                      << order << 
00106                      " not available for triangles");
00107 #endif
00108     }
00109 
00110   for (int i=0; i<q.length(); i++)
00111     {
00112       Array<Array<double> > qPerm;
00113       permute(multiplicity[i], q[i], qPerm);
00114       for (int j=0; j<multiplicity[i]; j++)
00115         {
00116           x.append(qPerm[j][0]);
00117           y.append(qPerm[j][1]);
00118           z.append(qPerm[j][2]);
00119           wgt.append(w[i]);
00120         }
00121     }
00122 }
00123 
00124 bool TetQuadrature::supportsOrder(int order)
00125 {
00126   if (order==1 || order==2 || order==4 || order==6) return true;
00127   return false;
00128 }
00129 
00130 void TetQuadrature::permute(int m, const Array<double>& q,
00131                             Array<Array<double> >& qPerm)
00132 {
00133   qPerm.resize(m);
00134   if (m==1)
00135     {
00136       qPerm[0] = tuple(q[0], q[0], q[0], q[0]);
00137     }
00138   else if (m==4)
00139     {
00140       qPerm[0] = tuple(q[0], q[1], q[1], q[1]);
00141       qPerm[1] = tuple(q[1], q[0], q[1], q[1]);
00142       qPerm[2] = tuple(q[1], q[1], q[0], q[1]);
00143       qPerm[3] = tuple(q[1], q[1], q[1], q[0]);
00144     }
00145   else if (m==12)
00146     {
00147       qPerm[0] = tuple(q[0], q[1], q[2], q[2]);
00148       qPerm[1] = tuple(q[0], q[2], q[1], q[2]);
00149       qPerm[2] = tuple(q[0], q[2], q[2], q[1]);
00150       qPerm[3] = tuple(q[1], q[0], q[2], q[2]);
00151       qPerm[4] = tuple(q[2], q[0], q[1], q[2]);
00152       qPerm[5] = tuple(q[2], q[0], q[2], q[1]);
00153       qPerm[6] = tuple(q[1], q[2], q[0], q[2]);
00154       qPerm[7] = tuple(q[2], q[1], q[0], q[2]);
00155       qPerm[8] = tuple(q[2], q[2], q[0], q[1]);
00156       qPerm[9] = tuple(q[1], q[2], q[2], q[0]);
00157       qPerm[10] = tuple(q[2], q[1], q[2], q[0]);
00158       qPerm[11] = tuple(q[2], q[2], q[1], q[0]);
00159     }
00160   else
00161     {
00162 #ifndef TRILINOS_7
00163       SUNDANCE_ERROR("invalid multiplicity " 
00164                      << m <<
00165                      " in TetQuadrature::permute()");
00166 #else
00167       SUNDANCE_ERROR7("invalid multiplicity " 
00168                      << m <<
00169                      " in TetQuadrature::permute()");
00170 #endif
00171     }
00172 }
00173 
00174 bool TetQuadrature::test(int p)
00175 {
00176   Array<double> w;
00177   Array<double> x;
00178   Array<double> y;
00179   Array<double> z;
00180 
00181   getPoints(p, w, x, y, z);
00182   bool pass = true;
00183   
00184   for (int a=0; a<=p; a++)
00185     {
00186       for (int b=0; b<p-a; b++)
00187         {
00188           int cMax = p - a - b;
00189           for (int c=0; c<=cMax; c++)
00190             {
00191               int dMax = p - a - b - c;
00192               for (int d=0; d<=dMax; d++)
00193                 {
00194                   double sum = 0.0;
00195                   for (int q=0; q<w.length(); q++)
00196                     {
00197                       sum += (1.0/6.0)*w[q] * pow(x[q], (double) a) * pow(y[q], (double) b)
00198                         * pow(z[q], (double) c)
00199                         * pow(1.0 - x[q] - y[q] - z[q], (double) d);
00200                     }
00201                   double err = fabs(sum - exact(a,b,c,d));
00202                   bool localPass = err < 1.0e-14;
00203                   pass = pass && localPass;
00204                   if (!localPass)
00205                     {
00206                       fprintf(stderr, "order=%d m (%d, %d, %d %d) q=%22.15g exact=%22.15g\n", p, a, b, c, d, sum, exact(a, b, c, d));
00207                       std::cerr << "error = " << err << std::endl;
00208                     }
00209                 }
00210             }
00211         }
00212     }
00213   return pass;
00214 }
00215 
00216 double TetQuadrature::exact(int a, int b, int c, int d)
00217 {
00218   return fact(a)*fact(b)*fact(c)*fact(d)/fact(a+b+c+d+3);
00219 }
00220 
00221 double TetQuadrature::fact(int x)
00222 {
00223   if (x==0) return 1.0;
00224   return x*fact(x-1);
00225 }
00226 

Site Contact