$extrastylesheet
quadrature_gm.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 #include "libmesh/quadrature_gm.h"
00020 
00021 namespace libMesh
00022 {
00023 
00024 // See also the file:
00025 // quadrature_gm_3D.C
00026 // for additional implementation.
00027 
00028 
00029 
00030 
00031 // Constructor
00032 QGrundmann_Moller::QGrundmann_Moller(const unsigned int d,
00033                                      const Order o) : QBase(d,o)
00034 {
00035 }
00036 
00037 
00038 // Destructor
00039 QGrundmann_Moller::~QGrundmann_Moller()
00040 {
00041 }
00042 
00043 
00044 
00045 
00046 
00047 void QGrundmann_Moller::gm_rule(unsigned int s, unsigned int dim)
00048 {
00049   // Only dim==2 or dim==3 is allowed
00050   libmesh_assert(dim==2 || dim==3);
00051 
00052   // A GM rule of index s can integrate polynomials of degree 2*s+1 exactly
00053   const unsigned int degree = 2*s+1;
00054 
00055   // The number of points for rule of index s is
00056   // (dim+1+s)! / (dim+1)! / s!
00057   // In 3D, this is = 1/24 * Pi_{i=1}^4 (s+i)
00058   // In 2D, this is =  1/6 * Pi_{i=1}^3 (s+i)
00059   const unsigned int n_pts = dim==2 ? (s+3)*(s+2)*(s+1) / 6 : (s+4)*(s+3)*(s+2)*(s+1) / 24;
00060   //libMesh::out << "n_pts=" << n_pts << std::endl;
00061 
00062   // Allocate space for points and weights
00063   _points.resize(n_pts);
00064   _weights.resize(n_pts);
00065 
00066   // (-1)^i -> This one flips sign at each iteration of the i-loop below.
00067   int one_pm=1;
00068 
00069   // Where we store all the integer point compositions/permutations
00070   std::vector<std::vector<unsigned int> > permutations;
00071 
00072   // Index into the vector where we should start adding the next round of points/weights
00073   std::size_t offset=0;
00074 
00075   // Implement the GM formula 4.1 on page 286 of the paper
00076   for (unsigned int i=0; i<=s; ++i)
00077     {
00078       // Get all the ordered compositions (and their permutations)
00079       // of |beta| = s-i into dim+1 parts
00080       compose_all(s-i, dim+1, permutations);
00081       //libMesh::out << "n. permutations=" << permutations.size() << std::endl;
00082 
00083       for (unsigned int p=0; p<permutations.size(); ++p)
00084         {
00085           // We use the first dim entries of each permutation to
00086           // construct an integration point.
00087           for (unsigned int j=0; j<dim; ++j)
00088             _points[offset+p](j) =
00089               static_cast<Real>(2.*permutations[p][j] + 1.) /
00090               static_cast<Real>(  degree + dim - 2.*i     );
00091         }
00092 
00093       // Compute the weight for this i, being careful to avoid overflow.
00094       // This technique is borrowed from Burkardt's code as well.
00095       // Use once for each of the points obtained from the permutations array.
00096       Real weight = one_pm;
00097 
00098       // This for loop needs to run for dim, degree, or dim+degree-i iterations,
00099       // whichever is largest.
00100       const unsigned int weight_loop_index =
00101         std::max(dim, std::max(degree, degree+dim-i))+1;
00102 
00103       for (unsigned int j=1; j<weight_loop_index; ++j)
00104         {
00105           if (j <= degree) // Accumulate (d+n-2i)^d term
00106             weight *= static_cast<Real>(degree+dim-2*i);
00107 
00108           if (j <= 2*s) // Accumulate 2^{-2s}
00109             weight *= 0.5;
00110 
00111           if (j <= i) // Accumulate (i!)^{-1}
00112             weight /= static_cast<Real>(j);
00113 
00114           if (j <= degree+dim-i) // Accumulate ( (d+n-i)! )^{-1}
00115             weight /= static_cast<Real>(j);
00116         }
00117 
00118       // This is the weight for each of the points computed previously
00119       for (unsigned int j=0; j<permutations.size(); ++j)
00120         _weights[offset+j] = weight;
00121 
00122       // Change sign for next iteration
00123       one_pm = -one_pm;
00124 
00125       // Update offset for the next set of points
00126       offset += permutations.size();
00127     }
00128 }
00129 
00130 
00131 
00132 
00133 // This routine for computing compositions and their permutations is
00134 // originall due to:
00135 //
00136 // Albert Nijenhuis, Herbert Wilf,
00137 // Combinatorial Algorithms for Computers and Calculators,
00138 // Second Edition,
00139 // Academic Press, 1978,
00140 // ISBN: 0-12-519260-6,
00141 // LC: QA164.N54.
00142 //
00143 // The routine is deceptively simple: I still don't understand exactly
00144 // why it works, but it does.
00145 void QGrundmann_Moller::compose_all(unsigned int s, // number to be compositioned
00146                                     unsigned int p, // # of partitions
00147                                     std::vector<std::vector<unsigned int> >& result)
00148 {
00149   // Clear out results remaining from previous calls
00150   result.clear();
00151 
00152   // Allocate storage for a workspace.  The workspace will periodically
00153   // be copied into the result container.
00154   std::vector<unsigned int> workspace(p);
00155 
00156   // The first result is always (s,0,...,0)
00157   workspace[0] = s;
00158   result.push_back(workspace);
00159 
00160   // the value of the first non-zero entry
00161   unsigned int head_value=s;
00162 
00163   // When head_index=-1, it refers to "off the front" of the array.  Therefore,
00164   // this needs to be a regular int rather than unsigned.  I initially tried to
00165   // do this with head_index unsigned and an else statement below, but then there
00166   // is the special case: (1,0,...,0) which does not work correctly.
00167   int head_index = -1;
00168 
00169   // At the end, all the entries will be in the final slot of workspace
00170   while (workspace.back() != s)
00171     {
00172       // Uncomment for debugging
00173       //libMesh::out << "previous head_value=" << head_value << " -> ";
00174 
00175       // If the previous head value is still larger than 1, reset the index
00176       // to "off the front" of the array
00177       if (head_value > 1)
00178         head_index = -1;
00179 
00180       // Either move the index onto the front of the array or on to
00181       // the next value.
00182       head_index++;
00183 
00184       // Get current value of the head entry
00185       head_value = workspace[head_index];
00186 
00187       // Uncomment for debugging
00188       //std::copy(workspace.begin(), workspace.end(), std::ostream_iterator<int>(libMesh::out, " "));
00189       //libMesh::out << ", head_index=" << head_index;
00190       //libMesh::out << ", head_value=" << head_value << " -> ";
00191 
00192       // Put a zero into the head_index of the array.  If head_index==0,
00193       // this will be overwritten in the next line with head_value-1.
00194       workspace[head_index] = 0;
00195 
00196       // The initial entry gets the current head value, minus 1.
00197       // If head_value > 1, the next loop iteration will start back
00198       // at workspace[0] again.
00199       libmesh_assert_greater (head_value, 0);
00200       workspace[0] = head_value - 1;
00201 
00202       // Increment the head+1 value
00203       workspace[head_index+1] += 1;
00204 
00205       // Save this composition in the results
00206       result.push_back(workspace);
00207 
00208       // Uncomment for debugging
00209       //std::copy(workspace.begin(), workspace.end(), std::ostream_iterator<int>(libMesh::out, " "));
00210       //libMesh::out<<"\n";
00211     }
00212 }
00213 
00214 } // namespace libMesh