$extrastylesheet
quadrature_monomial_2D.C
Go to the documentation of this file.
00001 // The libMesh Finite Element Library.
00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
00003 
00004 // This library is free software; you can redistribute it and/or
00005 // modify it under the terms of the GNU Lesser General Public
00006 // License as published by the Free Software Foundation; either
00007 // version 2.1 of the License, or (at your option) any later version.
00008 
00009 // This library is distributed in the hope that it will be useful,
00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00012 // Lesser General Public License for more details.
00013 
00014 // You should have received a copy of the GNU Lesser General Public
00015 // License along with this library; if not, write to the Free Software
00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00017 
00018 
00019 
00020 // Local includes
00021 #include "libmesh/quadrature_monomial.h"
00022 #include "libmesh/quadrature_gauss.h"
00023 
00024 namespace libMesh
00025 {
00026 
00027 
00028 void QMonomial::init_2D(const ElemType type_in,
00029                         unsigned int p)
00030 {
00031 
00032   switch (type_in)
00033     {
00034       //---------------------------------------------
00035       // Quadrilateral quadrature rules
00036     case QUAD4:
00037     case QUAD8:
00038     case QUAD9:
00039       {
00040         switch(_order + 2*p)
00041           {
00042           case SECOND:
00043             {
00044               // A degree=2 rule for the QUAD with 3 points.
00045               // A tensor product degree-2 Gauss would have 4 points.
00046               // This rule (or a variation on it) is probably available in
00047               //
00048               // A.H. Stroud, Approximate calculation of multiple integrals,
00049               // Prentice-Hall, Englewood Cliffs, N.J., 1971.
00050               //
00051               // though I have never actually seen a reference for it.
00052               // Luckily it's fairly easy to derive, which is what I've done
00053               // here [JWP].
00054               const Real
00055                 s=std::sqrt(1./3.),
00056                 t=std::sqrt(2./3.);
00057 
00058               const Real data[2][3] =
00059                 {
00060                   {0.0,  s,  2.0},
00061                   {  t, -s,  1.0}
00062                 };
00063 
00064               _points.resize(3);
00065               _weights.resize(3);
00066 
00067               wissmann_rule(data, 2);
00068 
00069               return;
00070             } // end case SECOND
00071 
00072 
00073 
00074             // For third-order, fall through to default case, use 2x2 Gauss product rule.
00075             // case THIRD:
00076             //   {
00077             //   }  // end case THIRD
00078 
00079           case FOURTH:
00080             {
00081               // A pair of degree=4 rules for the QUAD "C2" due to
00082               // Wissmann and Becker. These rules both have six points.
00083               // A tensor product degree-4 Gauss would have 9 points.
00084               //
00085               // J. W. Wissmann and T. Becker, Partially symmetric cubature
00086               // formulas for even degrees of exactness, SIAM J. Numer. Anal.  23
00087               // (1986), 676--685.
00088               const Real data[4][3] =
00089                 {
00090                   // First of 2 degree-4 rules given by Wissmann
00091                   {0.0000000000000000e+00,  0.0000000000000000e+00,  1.1428571428571428e+00},
00092                   {0.0000000000000000e+00,  9.6609178307929590e-01,  4.3956043956043956e-01},
00093                   {8.5191465330460049e-01,  4.5560372783619284e-01,  5.6607220700753210e-01},
00094                   {6.3091278897675402e-01, -7.3162995157313452e-01,  6.4271900178367668e-01}
00095                   //
00096                   // Second of 2 degree-4 rules given by Wissmann.  These both
00097                   // yield 4th-order accurate rules, I just chose the one that
00098                   // happened to contain the origin.
00099                   // {0.000000000000000, -0.356822089773090,  1.286412084888852},
00100                   // {0.000000000000000,  0.934172358962716,  0.491365692888926},
00101                   // {0.774596669241483,  0.390885162530071,  0.761883709085613},
00102                   // {0.774596669241483, -0.852765377881771,  0.349227402025498}
00103                 };
00104 
00105               _points.resize(6);
00106               _weights.resize(6);
00107 
00108               wissmann_rule(data, 4);
00109 
00110               return;
00111             } // end case FOURTH
00112 
00113 
00114 
00115 
00116           case FIFTH:
00117             {
00118               // A degree 5, 7-point rule due to Stroud.
00119               //
00120               // A.H. Stroud, Approximate calculation of multiple integrals,
00121               // Prentice-Hall, Englewood Cliffs, N.J., 1971.
00122               //
00123               // This rule is provably minimal in the number of points.
00124               // A tensor-product rule accurate for "bi-quintic" polynomials would have 9 points.
00125               const Real data[3][3] =
00126                 {
00127                   {                                  0.L,                                     0.L, static_cast<Real>(8.L  /  7.L)}, // 1
00128                   {                                  0.L, static_cast<Real>(std::sqrt(14.L/15.L)), static_cast<Real>(20.L / 63.L)}, // 2
00129                   {static_cast<Real>(std::sqrt(3.L/5.L)),   static_cast<Real>(std::sqrt(1.L/3.L)), static_cast<Real>(20.L / 36.L)}  // 4
00130                 };
00131 
00132               const unsigned int symmetry[3] = {
00133                 0, // Origin
00134                 7, // Central Symmetry
00135                 6  // Rectangular
00136               };
00137 
00138               _points.resize (7);
00139               _weights.resize(7);
00140 
00141               stroud_rule(data, symmetry, 3);
00142 
00143               return;
00144             } // end case FIFTH
00145 
00146 
00147 
00148 
00149           case SIXTH:
00150             {
00151               // A pair of degree=6 rules for the QUAD "C2" due to
00152               // Wissmann and Becker. These rules both have 10 points.
00153               // A tensor product degree-6 Gauss would have 16 points.
00154               //
00155               // J. W. Wissmann and T. Becker, Partially symmetric cubature
00156               // formulas for even degrees of exactness, SIAM J. Numer. Anal.  23
00157               // (1986), 676--685.
00158               const Real data[6][3] =
00159                 {
00160                   // First of 2 degree-6, 10 point rules given by Wissmann
00161                   // {0.000000000000000,  0.836405633697626,  0.455343245714174},
00162                   // {0.000000000000000, -0.357460165391307,  0.827395973202966},
00163                   // {0.888764014654765,  0.872101531193131,  0.144000884599645},
00164                   // {0.604857639464685,  0.305985162155427,  0.668259104262665},
00165                   // {0.955447506641064, -0.410270899466658,  0.225474004890679},
00166                   // {0.565459993438754, -0.872869311156879,  0.320896396788441}
00167                   //
00168                   // Second of 2 degree-6, 10 point rules given by Wissmann.
00169                   // Either of these will work, I just chose the one with points
00170                   // slightly further into the element interior.
00171                   {0.0000000000000000e+00,  8.6983337525005900e-01,  3.9275059096434794e-01},
00172                   {0.0000000000000000e+00, -4.7940635161211124e-01,  7.5476288124261053e-01},
00173                   {8.6374282634615388e-01,  8.0283751620765670e-01,  2.0616605058827902e-01},
00174                   {5.1869052139258234e-01,  2.6214366550805818e-01,  6.8999213848986375e-01},
00175                   {9.3397254497284950e-01, -3.6309658314806653e-01,  2.6051748873231697e-01},
00176                   {6.0897753601635630e-01, -8.9660863276245265e-01,  2.6956758608606100e-01}
00177                 };
00178 
00179               _points.resize(10);
00180               _weights.resize(10);
00181 
00182               wissmann_rule(data, 6);
00183 
00184               return;
00185             } // end case SIXTH
00186 
00187 
00188 
00189 
00190           case SEVENTH:
00191             {
00192               // A degree 7, 12-point rule due to Tyler, can be found in Stroud's book
00193               //
00194               // A.H. Stroud, Approximate calculation of multiple integrals,
00195               // Prentice-Hall, Englewood Cliffs, N.J., 1971.
00196               //
00197               // This rule is fully-symmetric and provably minimal in the number of points.
00198               // A tensor-product rule accurate for "bi-septic" polynomials would have 16 points.
00199               const Real
00200                 r  = std::sqrt(6.L/7.L),
00201                 s  = std::sqrt( (114.L - 3.L*std::sqrt(583.L)) / 287.L ),
00202                 t  = std::sqrt( (114.L + 3.L*std::sqrt(583.L)) / 287.L ),
00203                 B1 = 196.L / 810.L,
00204                 B2 = 4.L * (178981.L + 2769.L*std::sqrt(583.L)) / 1888920.L,
00205                 B3 = 4.L * (178981.L - 2769.L*std::sqrt(583.L)) / 1888920.L;
00206 
00207               const Real data[3][3] =
00208                 {
00209                   {r, 0.0, B1}, // 4
00210                   {s, 0.0, B2}, // 4
00211                   {t, 0.0, B3}  // 4
00212                 };
00213 
00214               const unsigned int symmetry[3] = {
00215                 3, // Full Symmetry, (x,0)
00216                 2, // Full Symmetry, (x,x)
00217                 2  // Full Symmetry, (x,x)
00218               };
00219 
00220               _points.resize (12);
00221               _weights.resize(12);
00222 
00223               stroud_rule(data, symmetry, 3);
00224 
00225               return;
00226             } // end case SEVENTH
00227 
00228 
00229 
00230 
00231           case EIGHTH:
00232             {
00233               // A pair of degree=8 rules for the QUAD "C2" due to
00234               // Wissmann and Becker. These rules both have 16 points.
00235               // A tensor product degree-6 Gauss would have 25 points.
00236               //
00237               // J. W. Wissmann and T. Becker, Partially symmetric cubature
00238               // formulas for even degrees of exactness, SIAM J. Numer. Anal.  23
00239               // (1986), 676--685.
00240               const Real data[10][3] =
00241                 {
00242                   // First of 2 degree-8, 16 point rules given by Wissmann
00243                   // {0.000000000000000,  0.000000000000000,  0.055364705621440},
00244                   // {0.000000000000000,  0.757629177660505,  0.404389368726076},
00245                   // {0.000000000000000, -0.236871842255702,  0.533546604952635},
00246                   // {0.000000000000000, -0.989717929044527,  0.117054188786739},
00247                   // {0.639091304900370,  0.950520955645667,  0.125614417613747},
00248                   // {0.937069076924990,  0.663882736885633,  0.136544584733588},
00249                   // {0.537083530541494,  0.304210681724104,  0.483408479211257},
00250                   // {0.887188506449625, -0.236496718536120,  0.252528506429544},
00251                   // {0.494698820670197, -0.698953476086564,  0.361262323882172},
00252                   // {0.897495818279768, -0.900390774211580,  0.085464254086247}
00253                   //
00254                   // Second of 2 degree-8, 16 point rules given by Wissmann.
00255                   // Either of these will work, I just chose the one with points
00256                   // further into the element interior.
00257                   {0.0000000000000000e+00,  6.5956013196034176e-01,  4.5027677630559029e-01},
00258                   {0.0000000000000000e+00, -9.4914292304312538e-01,  1.6657042677781274e-01},
00259                   {9.5250946607156228e-01,  7.6505181955768362e-01,  9.8869459933431422e-02},
00260                   {5.3232745407420624e-01,  9.3697598108841598e-01,  1.5369674714081197e-01},
00261                   {6.8473629795173504e-01,  3.3365671773574759e-01,  3.9668697607290278e-01},
00262                   {2.3314324080140552e-01, -7.9583272377396852e-02,  3.5201436794569501e-01},
00263                   {9.2768331930611748e-01, -2.7224008061253425e-01,  1.8958905457779799e-01},
00264                   {4.5312068740374942e-01, -6.1373535339802760e-01,  3.7510100114758727e-01},
00265                   {8.3750364042281223e-01, -8.8847765053597136e-01,  1.2561879164007201e-01}
00266                 };
00267 
00268               _points.resize(16);
00269               _weights.resize(16);
00270 
00271               wissmann_rule(data, /*10*/ 9);
00272 
00273               return;
00274             } // end case EIGHTH
00275 
00276 
00277 
00278 
00279           case NINTH:
00280             {
00281               // A degree 9, 17-point rule due to Moller.
00282               //
00283               // H.M. Moller,  Kubaturformeln mit minimaler Knotenzahl,
00284               // Numer. Math.  25 (1976), 185--200.
00285               //
00286               // This rule is provably minimal in the number of points.
00287               // A tensor-product rule accurate for "bi-ninth" degree polynomials would have 25 points.
00288               const Real data[5][3] =
00289                 {
00290                   {0.0000000000000000e+00, 0.0000000000000000e+00, 5.2674897119341563e-01}, // 1
00291                   {6.3068011973166885e-01, 9.6884996636197772e-01, 8.8879378170198706e-02}, // 4
00292                   {9.2796164595956966e-01, 7.5027709997890053e-01, 1.1209960212959648e-01}, // 4
00293                   {4.5333982113564719e-01, 5.2373582021442933e-01, 3.9828243926207009e-01}, // 4
00294                   {8.5261572933366230e-01, 7.6208328192617173e-02, 2.6905133763978080e-01}  // 4
00295                 };
00296 
00297               const unsigned int symmetry[5] = {
00298                 0, // Single point
00299                 4, // Rotational Invariant
00300                 4, // Rotational Invariant
00301                 4, // Rotational Invariant
00302                 4  // Rotational Invariant
00303               };
00304 
00305               _points.resize (17);
00306               _weights.resize(17);
00307 
00308               stroud_rule(data, symmetry, 5);
00309 
00310               return;
00311             } // end case NINTH
00312 
00313 
00314 
00315 
00316           case TENTH:
00317           case ELEVENTH:
00318             {
00319               // A degree 11, 24-point rule due to Cools and Haegemans.
00320               //
00321               // R. Cools and A. Haegemans, Another step forward in searching for
00322               // cubature formulae with a minimal number of knots for the square,
00323               // Computing 40 (1988), 139--146.
00324               //
00325               // P. Verlinden and R. Cools, The algebraic construction of a minimal
00326               // cubature formula of degree 11 for the square, Cubature Formulas
00327               // and their Applications (Russian) (Krasnoyarsk) (M.V. Noskov, ed.),
00328               // 1994, pp. 13--23.
00329               //
00330               // This rule is provably minimal in the number of points.
00331               // A tensor-product rule accurate for "bi-tenth" or "bi-eleventh" degree polynomials would have 36 points.
00332               const Real data[6][3] =
00333                 {
00334                   {6.9807610454956756e-01, 9.8263922354085547e-01, 4.8020763350723814e-02}, // 4
00335                   {9.3948638281673690e-01, 8.2577583590296393e-01, 6.6071329164550595e-02}, // 4
00336                   {9.5353952820153201e-01, 1.8858613871864195e-01, 9.7386777358668164e-02}, // 4
00337                   {3.1562343291525419e-01, 8.1252054830481310e-01, 2.1173634999894860e-01}, // 4
00338                   {7.1200191307533630e-01, 5.2532025036454776e-01, 2.2562606172886338e-01}, // 4
00339                   {4.2484724884866925e-01, 4.1658071912022368e-02, 3.5115871839824543e-01}  // 4
00340                 };
00341 
00342               const unsigned int symmetry[6] = {
00343                 4, // Rotational Invariant
00344                 4, // Rotational Invariant
00345                 4, // Rotational Invariant
00346                 4, // Rotational Invariant
00347                 4, // Rotational Invariant
00348                 4  // Rotational Invariant
00349               };
00350 
00351               _points.resize (24);
00352               _weights.resize(24);
00353 
00354               stroud_rule(data, symmetry, 6);
00355 
00356               return;
00357             } // end case TENTH,ELEVENTH
00358 
00359 
00360 
00361 
00362           case TWELFTH:
00363           case THIRTEENTH:
00364             {
00365               // A degree 13, 33-point rule due to Cools and Haegemans.
00366               //
00367               // R. Cools and A. Haegemans, Another step forward in searching for
00368               // cubature formulae with a minimal number of knots for the square,
00369               // Computing 40 (1988), 139--146.
00370               //
00371               // A tensor-product rule accurate for "bi-12" or "bi-13" degree polynomials would have 49 points.
00372               const Real data[9][3] =
00373                 {
00374                   {0.0000000000000000e+00, 0.0000000000000000e+00, 3.0038211543122536e-01}, // 1
00375                   {9.8348668243987226e-01, 7.7880971155441942e-01, 2.9991838864499131e-02}, // 4
00376                   {8.5955600564163892e-01, 9.5729769978630736e-01, 3.8174421317083669e-02}, // 4
00377                   {9.5892517028753485e-01, 1.3818345986246535e-01, 6.0424923817749980e-02}, // 4
00378                   {3.9073621612946100e-01, 9.4132722587292523e-01, 7.7492738533105339e-02}, // 4
00379                   {8.5007667369974857e-01, 4.7580862521827590e-01, 1.1884466730059560e-01}, // 4
00380                   {6.4782163718701073e-01, 7.5580535657208143e-01, 1.2976355037000271e-01}, // 4
00381                   {7.0741508996444936e-02, 6.9625007849174941e-01, 2.1334158145718938e-01}, // 4
00382                   {4.0930456169403884e-01, 3.4271655604040678e-01, 2.5687074948196783e-01}  // 4
00383                 };
00384 
00385               const unsigned int symmetry[9] = {
00386                 0, // Single point
00387                 4, // Rotational Invariant
00388                 4, // Rotational Invariant
00389                 4, // Rotational Invariant
00390                 4, // Rotational Invariant
00391                 4, // Rotational Invariant
00392                 4, // Rotational Invariant
00393                 4, // Rotational Invariant
00394                 4  // Rotational Invariant
00395               };
00396 
00397               _points.resize (33);
00398               _weights.resize(33);
00399 
00400               stroud_rule(data, symmetry, 9);
00401 
00402               return;
00403             } // end case TWELFTH,THIRTEENTH
00404 
00405 
00406 
00407 
00408           case FOURTEENTH:
00409           case FIFTEENTH:
00410             {
00411               // A degree-15, 48 point rule originally due to Rabinowitz and Richter,
00412               // can be found in Cools' 1971 book.
00413               //
00414               // A.H. Stroud, Approximate calculation of multiple integrals,
00415               // Prentice-Hall, Englewood Cliffs, N.J., 1971.
00416               //
00417               // The product Gauss rule for this order has 8^2=64 points.
00418               const Real data[9][3] =
00419                 {
00420                   {9.915377816777667e-01L, 0.0000000000000000e+00,  3.01245207981210e-02L}, // 4
00421                   {8.020163879230440e-01L, 0.0000000000000000e+00,  8.71146840209092e-02L}, // 4
00422                   {5.648674875232742e-01L, 0.0000000000000000e+00, 1.250080294351494e-01L}, // 4
00423                   {9.354392392539896e-01L, 0.0000000000000000e+00,  2.67651407861666e-02L}, // 4
00424                   {7.624563338825799e-01L, 0.0000000000000000e+00,  9.59651863624437e-02L}, // 4
00425                   {2.156164241427213e-01L, 0.0000000000000000e+00, 1.750832998343375e-01L}, // 4
00426                   {9.769662659711761e-01L, 6.684480048977932e-01L,  2.83136372033274e-02L}, // 4
00427                   {8.937128379503403e-01L, 3.735205277617582e-01L,  8.66414716025093e-02L}, // 4
00428                   {6.122485619312083e-01L, 4.078983303613935e-01L, 1.150144605755996e-01L}  // 4
00429                 };
00430 
00431               const unsigned int symmetry[9] = {
00432                 3, // Full Symmetry, (x,0)
00433                 3, // Full Symmetry, (x,0)
00434                 3, // Full Symmetry, (x,0)
00435                 2, // Full Symmetry, (x,x)
00436                 2, // Full Symmetry, (x,x)
00437                 2, // Full Symmetry, (x,x)
00438                 1, // Full Symmetry, (x,y)
00439                 1, // Full Symmetry, (x,y)
00440                 1, // Full Symmetry, (x,y)
00441               };
00442 
00443               _points.resize (48);
00444               _weights.resize(48);
00445 
00446               stroud_rule(data, symmetry, 9);
00447 
00448               return;
00449             } //   case FOURTEENTH, FIFTEENTH:
00450 
00451 
00452 
00453 
00454           case SIXTEENTH:
00455           case SEVENTEENTH:
00456             {
00457               // A degree 17, 60-point rule due to Cools and Haegemans.
00458               //
00459               // R. Cools and A. Haegemans, Another step forward in searching for
00460               // cubature formulae with a minimal number of knots for the square,
00461               // Computing 40 (1988), 139--146.
00462               //
00463               // A tensor-product rule accurate for "bi-14" or "bi-15" degree polynomials would have 64 points.
00464               // A tensor-product rule accurate for "bi-16" or "bi-17" degree polynomials would have 81 points.
00465               const Real data[10][3] =
00466                 {
00467                   {9.8935307451260049e-01, 0.0000000000000000e+00, 2.0614915919990959e-02}, // 4
00468                   {3.7628520715797329e-01, 0.0000000000000000e+00, 1.2802571617990983e-01}, // 4
00469                   {9.7884827926223311e-01, 0.0000000000000000e+00, 5.5117395340318905e-03}, // 4
00470                   {8.8579472916411612e-01, 0.0000000000000000e+00, 3.9207712457141880e-02}, // 4
00471                   {1.7175612383834817e-01, 0.0000000000000000e+00, 7.6396945079863302e-02}, // 4
00472                   {5.9049927380600241e-01, 3.1950503663457394e-01, 1.4151372994997245e-01}, // 8
00473                   {7.9907913191686325e-01, 5.9797245192945738e-01, 8.3903279363797602e-02}, // 8
00474                   {8.0374396295874471e-01, 5.8344481776550529e-02, 6.0394163649684546e-02}, // 8
00475                   {9.3650627612749478e-01, 3.4738631616620267e-01, 5.7387752969212695e-02}, // 8
00476                   {9.8132117980545229e-01, 7.0600028779864611e-01, 2.1922559481863763e-02}, // 8
00477                 };
00478 
00479               const unsigned int symmetry[10] = {
00480                 3, // Fully symmetric (x,0)
00481                 3, // Fully symmetric (x,0)
00482                 2, // Fully symmetric (x,x)
00483                 2, // Fully symmetric (x,x)
00484                 2, // Fully symmetric (x,x)
00485                 1, // Fully symmetric (x,y)
00486                 1, // Fully symmetric (x,y)
00487                 1, // Fully symmetric (x,y)
00488                 1, // Fully symmetric (x,y)
00489                 1  // Fully symmetric (x,y)
00490               };
00491 
00492               _points.resize (60);
00493               _weights.resize(60);
00494 
00495               stroud_rule(data, symmetry, 10);
00496 
00497               return;
00498             } // end case FOURTEENTH through SEVENTEENTH
00499 
00500 
00501 
00502             // By default: construct and use a Gauss quadrature rule
00503           default:
00504             {
00505               // Break out and fall down into the default: case for the
00506               // outer switch statement.
00507               break;
00508             }
00509 
00510           } // end switch(_order + 2*p)
00511       } // end case QUAD4/8/9
00512 
00513 
00514       // By default: construct and use a Gauss quadrature rule
00515     default:
00516       {
00517         QGauss gauss_rule(2, _order);
00518         gauss_rule.init(type_in, p);
00519 
00520         // Swap points and weights with the about-to-be destroyed rule.
00521         _points.swap (gauss_rule.get_points() );
00522         _weights.swap(gauss_rule.get_weights());
00523 
00524         return;
00525       }
00526     } // end switch (type_in)
00527 }
00528 
00529 } // namespace libMesh