quadrature_gauss.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_gauss.h" 00020 00021 namespace libMesh 00022 { 00023 00024 // See the files: 00025 00026 // quadrature_gauss_1D.C 00027 // quadrature_gauss_2D.C 00028 // quadrature_gauss_3D.C 00029 00030 // for implementation of specific element types. 00031 00032 00033 void QGauss::keast_rule(const Real rule_data[][4], 00034 const unsigned int n_pts) 00035 { 00036 // Like the Dunavant rule, the input data should have 4 columns. These columns 00037 // have the following format and implied permutations (w=weight). 00038 // {a, 0, 0, w} = 1-permutation (a,a,a) 00039 // {a, b, 0, w} = 4-permutation (a,b,b), (b,a,b), (b,b,a), (b,b,b) 00040 // {a, 0, b, w} = 6-permutation (a,a,b), (a,b,b), (b,b,a), (b,a,b), (b,a,a), (a,b,a) 00041 // {a, b, c, w} = 12-permutation (a,a,b), (a,a,c), (b,a,a), (c,a,a), (a,b,a), (a,c,a) 00042 // (a,b,c), (a,c,b), (b,a,c), (b,c,a), (c,a,b), (c,b,a) 00043 00044 // Always insert into the points & weights vector relative to the offset 00045 unsigned int offset=0; 00046 00047 00048 for (unsigned int p=0; p<n_pts; ++p) 00049 { 00050 00051 // There must always be a non-zero entry to start the row 00052 libmesh_assert_not_equal_to (rule_data[p][0], static_cast<Real>(0.0)); 00053 00054 // A zero weight may imply you did not set up the raw data correctly 00055 libmesh_assert_not_equal_to (rule_data[p][3], static_cast<Real>(0.0)); 00056 00057 // What kind of point is this? 00058 // One non-zero entry in first 3 cols ? 1-perm (centroid) point = 1 00059 // Two non-zero entries in first 3 cols ? 3-perm point = 3 00060 // Three non-zero entries ? 6-perm point = 6 00061 unsigned int pointtype=1; 00062 00063 if (rule_data[p][1] != static_cast<Real>(0.0)) 00064 { 00065 if (rule_data[p][2] != static_cast<Real>(0.0)) 00066 pointtype = 12; 00067 else 00068 pointtype = 4; 00069 } 00070 else 00071 { 00072 // The second entry is zero. What about the third? 00073 if (rule_data[p][2] != static_cast<Real>(0.0)) 00074 pointtype = 6; 00075 } 00076 00077 00078 switch (pointtype) 00079 { 00080 case 1: 00081 { 00082 // Be sure we have enough space to insert this point 00083 libmesh_assert_less (offset + 0, _points.size()); 00084 00085 const Real a = rule_data[p][0]; 00086 00087 // The point has only a single permutation (the centroid!) 00088 _points[offset + 0] = Point(a,a,a); 00089 00090 // The weight is always the last entry in the row. 00091 _weights[offset + 0] = rule_data[p][3]; 00092 00093 offset += pointtype; 00094 break; 00095 } 00096 00097 case 4: 00098 { 00099 // Be sure we have enough space to insert these points 00100 libmesh_assert_less (offset + 3, _points.size()); 00101 00102 const Real a = rule_data[p][0]; 00103 const Real b = rule_data[p][1]; 00104 const Real wt = rule_data[p][3]; 00105 00106 // Here it's understood the second entry is to be used twice, and 00107 // thus there are three possible permutations. 00108 _points[offset + 0] = Point(a,b,b); 00109 _points[offset + 1] = Point(b,a,b); 00110 _points[offset + 2] = Point(b,b,a); 00111 _points[offset + 3] = Point(b,b,b); 00112 00113 for (unsigned int j=0; j<pointtype; ++j) 00114 _weights[offset + j] = wt; 00115 00116 offset += pointtype; 00117 break; 00118 } 00119 00120 case 6: 00121 { 00122 // Be sure we have enough space to insert these points 00123 libmesh_assert_less (offset + 5, _points.size()); 00124 00125 const Real a = rule_data[p][0]; 00126 const Real b = rule_data[p][2]; 00127 const Real wt = rule_data[p][3]; 00128 00129 // Three individual entries with six permutations. 00130 _points[offset + 0] = Point(a,a,b); 00131 _points[offset + 1] = Point(a,b,b); 00132 _points[offset + 2] = Point(b,b,a); 00133 _points[offset + 3] = Point(b,a,b); 00134 _points[offset + 4] = Point(b,a,a); 00135 _points[offset + 5] = Point(a,b,a); 00136 00137 for (unsigned int j=0; j<pointtype; ++j) 00138 _weights[offset + j] = wt; 00139 00140 offset += pointtype; 00141 break; 00142 } 00143 00144 00145 case 12: 00146 { 00147 // Be sure we have enough space to insert these points 00148 libmesh_assert_less (offset + 11, _points.size()); 00149 00150 const Real a = rule_data[p][0]; 00151 const Real b = rule_data[p][1]; 00152 const Real c = rule_data[p][2]; 00153 const Real wt = rule_data[p][3]; 00154 00155 // Three individual entries with six permutations. 00156 _points[offset + 0] = Point(a,a,b); _points[offset + 6] = Point(a,b,c); 00157 _points[offset + 1] = Point(a,a,c); _points[offset + 7] = Point(a,c,b); 00158 _points[offset + 2] = Point(b,a,a); _points[offset + 8] = Point(b,a,c); 00159 _points[offset + 3] = Point(c,a,a); _points[offset + 9] = Point(b,c,a); 00160 _points[offset + 4] = Point(a,b,a); _points[offset + 10] = Point(c,a,b); 00161 _points[offset + 5] = Point(a,c,a); _points[offset + 11] = Point(c,b,a); 00162 00163 for (unsigned int j=0; j<pointtype; ++j) 00164 _weights[offset + j] = wt; 00165 00166 offset += pointtype; 00167 break; 00168 } 00169 00170 default: 00171 { 00172 libMesh::err << "Don't know what to do with this many permutation points!" << std::endl; 00173 libmesh_error(); 00174 } 00175 } 00176 00177 } 00178 00179 } 00180 00181 00182 // A number of different rules for triangles can be described by 00183 // permutations of the following types of points: 00184 // I: "1"-permutation, (1/3,1/3) (single point only) 00185 // II: 3-permutation, (a,a,1-2a) 00186 // III: 6-permutation, (a,b,1-a-b) 00187 // The weights for a given set of permutations are all the same. 00188 void QGauss::dunavant_rule2(const Real* wts, 00189 const Real* a, 00190 const Real* b, 00191 const unsigned int* permutation_ids, 00192 unsigned int n_wts) 00193 { 00194 // Figure out how many total points by summing up the entries 00195 // in the permutation_ids array, and resize the _points and _weights 00196 // vectors appropriately. 00197 unsigned int total_pts = 0; 00198 for (unsigned int p=0; p<n_wts; ++p) 00199 total_pts += permutation_ids[p]; 00200 00201 // Resize point and weight vectors appropriately. 00202 _points.resize(total_pts); 00203 _weights.resize(total_pts); 00204 00205 // Always insert into the points & weights vector relative to the offset 00206 unsigned int offset=0; 00207 00208 for (unsigned int p=0; p<n_wts; ++p) 00209 { 00210 switch (permutation_ids[p]) 00211 { 00212 case 1: 00213 { 00214 // The point has only a single permutation (the centroid!) 00215 // So we don't even need to look in the a or b arrays. 00216 _points [offset + 0] = Point(1.0L/3.0L, 1.0L/3.0L); 00217 _weights[offset + 0] = wts[p]; 00218 00219 offset += 1; 00220 break; 00221 } 00222 00223 00224 case 3: 00225 { 00226 // For this type of rule, don't need to look in the b array. 00227 _points[offset + 0] = Point(a[p], a[p]); // (a,a) 00228 _points[offset + 1] = Point(a[p], 1.L-2.L*a[p]); // (a,1-2a) 00229 _points[offset + 2] = Point(1.L-2.L*a[p], a[p]); // (1-2a,a) 00230 00231 for (unsigned int j=0; j<3; ++j) 00232 _weights[offset + j] = wts[p]; 00233 00234 offset += 3; 00235 break; 00236 } 00237 00238 case 6: 00239 { 00240 // This type of point uses all 3 arrays... 00241 _points[offset + 0] = Point(a[p], b[p]); 00242 _points[offset + 1] = Point(b[p], a[p]); 00243 _points[offset + 2] = Point(a[p], 1.L-a[p]-b[p]); 00244 _points[offset + 3] = Point(1.L-a[p]-b[p], a[p]); 00245 _points[offset + 4] = Point(b[p], 1.L-a[p]-b[p]); 00246 _points[offset + 5] = Point(1.L-a[p]-b[p], b[p]); 00247 00248 for (unsigned int j=0; j<6; ++j) 00249 _weights[offset + j] = wts[p]; 00250 00251 offset += 6; 00252 break; 00253 } 00254 00255 default: 00256 { 00257 libMesh::err << "Unknown permutation id: " << permutation_ids[p] << "!" << std::endl; 00258 libmesh_error(); 00259 } 00260 } 00261 } 00262 00263 } 00264 00265 00266 void QGauss::dunavant_rule(const Real rule_data[][4], 00267 const unsigned int n_pts) 00268 { 00269 // The input data array has 4 columns. The first 3 are the permutation points. 00270 // The last column is the weights for a given set of permutation points. A zero 00271 // in two of the first 3 columns implies the point is a 1-permutation (centroid). 00272 // A zero in one of the first 3 columns implies the point is a 3-permutation. 00273 // Otherwise each point is assumed to be a 6-permutation. 00274 00275 // Always insert into the points & weights vector relative to the offset 00276 unsigned int offset=0; 00277 00278 00279 for (unsigned int p=0; p<n_pts; ++p) 00280 { 00281 00282 // There must always be a non-zero entry to start the row 00283 libmesh_assert_not_equal_to ( rule_data[p][0], static_cast<Real>(0.0) ); 00284 00285 // A zero weight may imply you did not set up the raw data correctly 00286 libmesh_assert_not_equal_to ( rule_data[p][3], static_cast<Real>(0.0) ); 00287 00288 // What kind of point is this? 00289 // One non-zero entry in first 3 cols ? 1-perm (centroid) point = 1 00290 // Two non-zero entries in first 3 cols ? 3-perm point = 3 00291 // Three non-zero entries ? 6-perm point = 6 00292 unsigned int pointtype=1; 00293 00294 if (rule_data[p][1] != static_cast<Real>(0.0)) 00295 { 00296 if (rule_data[p][2] != static_cast<Real>(0.0)) 00297 pointtype = 6; 00298 else 00299 pointtype = 3; 00300 } 00301 00302 switch (pointtype) 00303 { 00304 case 1: 00305 { 00306 // Be sure we have enough space to insert this point 00307 libmesh_assert_less (offset + 0, _points.size()); 00308 00309 // The point has only a single permutation (the centroid!) 00310 _points[offset + 0] = Point(rule_data[p][0], rule_data[p][0]); 00311 00312 // The weight is always the last entry in the row. 00313 _weights[offset + 0] = rule_data[p][3]; 00314 00315 offset += 1; 00316 break; 00317 } 00318 00319 case 3: 00320 { 00321 // Be sure we have enough space to insert these points 00322 libmesh_assert_less (offset + 2, _points.size()); 00323 00324 // Here it's understood the second entry is to be used twice, and 00325 // thus there are three possible permutations. 00326 _points[offset + 0] = Point(rule_data[p][0], rule_data[p][1]); 00327 _points[offset + 1] = Point(rule_data[p][1], rule_data[p][0]); 00328 _points[offset + 2] = Point(rule_data[p][1], rule_data[p][1]); 00329 00330 for (unsigned int j=0; j<3; ++j) 00331 _weights[offset + j] = rule_data[p][3]; 00332 00333 offset += 3; 00334 break; 00335 } 00336 00337 case 6: 00338 { 00339 // Be sure we have enough space to insert these points 00340 libmesh_assert_less (offset + 5, _points.size()); 00341 00342 // Three individual entries with six permutations. 00343 _points[offset + 0] = Point(rule_data[p][0], rule_data[p][1]); 00344 _points[offset + 1] = Point(rule_data[p][0], rule_data[p][2]); 00345 _points[offset + 2] = Point(rule_data[p][1], rule_data[p][0]); 00346 _points[offset + 3] = Point(rule_data[p][1], rule_data[p][2]); 00347 _points[offset + 4] = Point(rule_data[p][2], rule_data[p][0]); 00348 _points[offset + 5] = Point(rule_data[p][2], rule_data[p][1]); 00349 00350 for (unsigned int j=0; j<6; ++j) 00351 _weights[offset + j] = rule_data[p][3]; 00352 00353 offset += 6; 00354 break; 00355 } 00356 00357 default: 00358 { 00359 libMesh::err << "Don't know what to do with this many permutation points!" << std::endl; 00360 libmesh_error(); 00361 } 00362 } 00363 } 00364 } 00365 00366 } // namespace libMesh 00367 00368
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:48 UTC
Hosted By: