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