$extrastylesheet
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_gauss.h" 00022 #include "libmesh/quadrature_conical.h" 00023 #include "libmesh/quadrature_gm.h" 00024 00025 namespace libMesh 00026 { 00027 00028 void QGauss::init_3D(const ElemType type_in, 00029 unsigned int p) 00030 { 00031 #if LIBMESH_DIM == 3 00032 00033 //----------------------------------------------------------------------- 00034 // 3D quadrature rules 00035 switch (type_in) 00036 { 00037 //--------------------------------------------- 00038 // Hex quadrature rules 00039 case HEX8: 00040 case HEX20: 00041 case HEX27: 00042 { 00043 // We compute the 3D quadrature rule as a tensor 00044 // product of the 1D quadrature rule. 00045 QGauss q1D(1,_order); 00046 q1D.init(EDGE2,p); 00047 00048 tensor_product_hex( q1D ); 00049 00050 return; 00051 } 00052 00053 00054 00055 //--------------------------------------------- 00056 // Tetrahedral quadrature rules 00057 case TET4: 00058 case TET10: 00059 { 00060 switch(_order + 2*p) 00061 { 00062 // Taken from pg. 222 of "The finite element method," vol. 1 00063 // ed. 5 by Zienkiewicz & Taylor 00064 case CONSTANT: 00065 case FIRST: 00066 { 00067 // Exact for linears 00068 _points.resize(1); 00069 _weights.resize(1); 00070 00071 00072 _points[0](0) = .25; 00073 _points[0](1) = .25; 00074 _points[0](2) = .25; 00075 00076 _weights[0] = .1666666666666666666666666666666666666666666667L; 00077 00078 return; 00079 } 00080 case SECOND: 00081 { 00082 // Exact for quadratics 00083 _points.resize(4); 00084 _weights.resize(4); 00085 00086 00087 const Real a = .585410196624969; 00088 const Real b = .138196601125011; 00089 00090 _points[0](0) = a; 00091 _points[0](1) = b; 00092 _points[0](2) = b; 00093 00094 _points[1](0) = b; 00095 _points[1](1) = a; 00096 _points[1](2) = b; 00097 00098 _points[2](0) = b; 00099 _points[2](1) = b; 00100 _points[2](2) = a; 00101 00102 _points[3](0) = b; 00103 _points[3](1) = b; 00104 _points[3](2) = b; 00105 00106 00107 00108 _weights[0] = .0416666666666666666666666666666666666666666667; 00109 _weights[1] = _weights[0]; 00110 _weights[2] = _weights[0]; 00111 _weights[3] = _weights[0]; 00112 00113 return; 00114 } 00115 00116 00117 00118 // Can be found in the class notes 00119 // http://www.cs.rpi.edu/~flaherje/FEM/fem6.ps 00120 // by Flaherty. 00121 // 00122 // Caution: this rule has a negative weight and may be 00123 // unsuitable for some problems. 00124 // Exact for cubics. 00125 // 00126 // Note: Keast (see ref. elsewhere in this file) also gives 00127 // a third-order rule with positive weights, but it contains points 00128 // on the ref. elt. boundary, making it less suitable for FEM calculations. 00129 case THIRD: 00130 { 00131 if (allow_rules_with_negative_weights) 00132 { 00133 _points.resize(5); 00134 _weights.resize(5); 00135 00136 00137 _points[0](0) = .25; 00138 _points[0](1) = .25; 00139 _points[0](2) = .25; 00140 00141 _points[1](0) = .5; 00142 _points[1](1) = .16666666666666666666666666666666666666666667; 00143 _points[1](2) = .16666666666666666666666666666666666666666667; 00144 00145 _points[2](0) = .16666666666666666666666666666666666666666667; 00146 _points[2](1) = .5; 00147 _points[2](2) = .16666666666666666666666666666666666666666667; 00148 00149 _points[3](0) = .16666666666666666666666666666666666666666667; 00150 _points[3](1) = .16666666666666666666666666666666666666666667; 00151 _points[3](2) = .5; 00152 00153 _points[4](0) = .16666666666666666666666666666666666666666667; 00154 _points[4](1) = .16666666666666666666666666666666666666666667; 00155 _points[4](2) = .16666666666666666666666666666666666666666667; 00156 00157 00158 _weights[0] = -.133333333333333333333333333333333333333333333; 00159 _weights[1] = .075; 00160 _weights[2] = _weights[1]; 00161 _weights[3] = _weights[1]; 00162 _weights[4] = _weights[1]; 00163 00164 return; 00165 } // end if (allow_rules_with_negative_weights) 00166 else 00167 { 00168 // If a rule with postive weights is required, the 2x2x2 conical 00169 // product rule is third-order accurate and has less points than 00170 // the next-available positive-weight rule at FIFTH order. 00171 QConical conical_rule(3, _order); 00172 conical_rule.init(type_in, p); 00173 00174 // Swap points and weights with the about-to-be destroyed rule. 00175 _points.swap (conical_rule.get_points() ); 00176 _weights.swap(conical_rule.get_weights()); 00177 00178 return; 00179 } 00180 // Note: if !allow_rules_with_negative_weights, fall through to next case. 00181 } 00182 00183 00184 00185 // Originally a Keast rule, 00186 // Patrick Keast, 00187 // Moderate Degree Tetrahedral Quadrature Formulas, 00188 // Computer Methods in Applied Mechanics and Engineering, 00189 // Volume 55, Number 3, May 1986, pages 339-348. 00190 // 00191 // Can also be found the class notes 00192 // http://www.cs.rpi.edu/~flaherje/FEM/fem6.ps 00193 // by Flaherty. 00194 // 00195 // Caution: this rule has a negative weight and may be 00196 // unsuitable for some problems. 00197 case FOURTH: 00198 { 00199 if (allow_rules_with_negative_weights) 00200 { 00201 _points.resize(11); 00202 _weights.resize(11); 00203 00204 // The raw data for the quadrature rule. 00205 const Real rule_data[3][4] = { 00206 {0.250000000000000000e+00, 0., 0., -0.131555555555555556e-01}, // 1 00207 {0.785714285714285714e+00, 0.714285714285714285e-01, 0., 0.762222222222222222e-02}, // 4 00208 {0.399403576166799219e+00, 0., 0.100596423833200785e+00, 0.248888888888888889e-01} // 6 00209 }; 00210 00211 00212 // Now call the keast routine to generate _points and _weights 00213 keast_rule(rule_data, 3); 00214 00215 return; 00216 } // end if (allow_rules_with_negative_weights) 00217 // Note: if !allow_rules_with_negative_weights, fall through to next case. 00218 } 00219 00220 00221 00222 00223 // Walkington's fifth-order 14-point rule from 00224 // "Quadrature on Simplices of Arbitrary Dimension" 00225 // 00226 // We originally had a Keast rule here, but this rule had 00227 // more points than an equivalent rule by Walkington and 00228 // also contained points on the boundary of the ref. elt, 00229 // making it less suitable for FEM calculations. 00230 case FIFTH: 00231 { 00232 _points.resize(14); 00233 _weights.resize(14); 00234 00235 // permutations of these points and suitably-modified versions of 00236 // these points are the quadrature point locations 00237 const Real a[3] = {0.31088591926330060980, // a1 from the paper 00238 0.092735250310891226402, // a2 from the paper 00239 0.045503704125649649492}; // a3 from the paper 00240 00241 // weights. a[] and wt[] are the only floating-point inputs required 00242 // for this rule. 00243 const Real wt[3] = {0.018781320953002641800, // w1 from the paper 00244 0.012248840519393658257, // w2 from the paper 00245 0.0070910034628469110730}; // w3 from the paper 00246 00247 // The first two sets of 4 points are formed in a similar manner 00248 for (unsigned int i=0; i<2; ++i) 00249 { 00250 // Where we will insert values into _points and _weights 00251 const unsigned int offset=4*i; 00252 00253 // Stuff points and weights values into their arrays 00254 const Real b = 1. - 3.*a[i]; 00255 00256 // Here are the permutations. Order of these is not important, 00257 // all have the same weight 00258 _points[offset + 0] = Point(a[i], a[i], a[i]); 00259 _points[offset + 1] = Point(a[i], b, a[i]); 00260 _points[offset + 2] = Point( b, a[i], a[i]); 00261 _points[offset + 3] = Point(a[i], a[i], b); 00262 00263 // These 4 points all have the same weights 00264 for (unsigned int j=0; j<4; ++j) 00265 _weights[offset + j] = wt[i]; 00266 } // end for 00267 00268 00269 { 00270 // The third set contains 6 points and is formed a little differently 00271 const unsigned int offset = 8; 00272 const Real b = 0.5*(1. - 2.*a[2]); 00273 00274 // Here are the permutations. Order of these is not important, 00275 // all have the same weight 00276 _points[offset + 0] = Point(b , b, a[2]); 00277 _points[offset + 1] = Point(b , a[2], a[2]); 00278 _points[offset + 2] = Point(a[2], a[2], b); 00279 _points[offset + 3] = Point(a[2], b, a[2]); 00280 _points[offset + 4] = Point( b, a[2], b); 00281 _points[offset + 5] = Point(a[2], b, b); 00282 00283 // These 6 points all have the same weights 00284 for (unsigned int j=0; j<6; ++j) 00285 _weights[offset + j] = wt[2]; 00286 } 00287 00288 00289 // Original rule by Keast, unsuitable because it has points on the 00290 // reference element boundary! 00291 // _points.resize(15); 00292 // _weights.resize(15); 00293 00294 // _points[0](0) = 0.25; 00295 // _points[0](1) = 0.25; 00296 // _points[0](2) = 0.25; 00297 00298 // { 00299 // const Real a = 0.; 00300 // const Real b = 0.333333333333333333333333333333333333333; 00301 00302 // _points[1](0) = a; 00303 // _points[1](1) = b; 00304 // _points[1](2) = b; 00305 00306 // _points[2](0) = b; 00307 // _points[2](1) = a; 00308 // _points[2](2) = b; 00309 00310 // _points[3](0) = b; 00311 // _points[3](1) = b; 00312 // _points[3](2) = a; 00313 00314 // _points[4](0) = b; 00315 // _points[4](1) = b; 00316 // _points[4](2) = b; 00317 // } 00318 // { 00319 // const Real a = 0.7272727272727272727272727272727272727272727272727272727; 00320 // const Real b = 0.0909090909090909090909090909090909090909090909090909091; 00321 00322 // _points[5](0) = a; 00323 // _points[5](1) = b; 00324 // _points[5](2) = b; 00325 00326 // _points[6](0) = b; 00327 // _points[6](1) = a; 00328 // _points[6](2) = b; 00329 00330 // _points[7](0) = b; 00331 // _points[7](1) = b; 00332 // _points[7](2) = a; 00333 00334 // _points[8](0) = b; 00335 // _points[8](1) = b; 00336 // _points[8](2) = b; 00337 // } 00338 // { 00339 // const Real a = 0.066550153573664; 00340 // const Real b = 0.433449846426336; 00341 00342 // _points[9](0) = b; 00343 // _points[9](1) = a; 00344 // _points[9](2) = a; 00345 00346 // _points[10](0) = a; 00347 // _points[10](1) = a; 00348 // _points[10](2) = b; 00349 00350 // _points[11](0) = a; 00351 // _points[11](1) = b; 00352 // _points[11](2) = b; 00353 00354 // _points[12](0) = b; 00355 // _points[12](1) = a; 00356 // _points[12](2) = b; 00357 00358 // _points[13](0) = b; 00359 // _points[13](1) = b; 00360 // _points[13](2) = a; 00361 00362 // _points[14](0) = a; 00363 // _points[14](1) = b; 00364 // _points[14](2) = a; 00365 // } 00366 00367 // _weights[0] = 0.030283678097089; 00368 // _weights[1] = 0.006026785714286; 00369 // _weights[2] = _weights[1]; 00370 // _weights[3] = _weights[1]; 00371 // _weights[4] = _weights[1]; 00372 // _weights[5] = 0.011645249086029; 00373 // _weights[6] = _weights[5]; 00374 // _weights[7] = _weights[5]; 00375 // _weights[8] = _weights[5]; 00376 // _weights[9] = 0.010949141561386; 00377 // _weights[10] = _weights[9]; 00378 // _weights[11] = _weights[9]; 00379 // _weights[12] = _weights[9]; 00380 // _weights[13] = _weights[9]; 00381 // _weights[14] = _weights[9]; 00382 00383 return; 00384 } 00385 00386 00387 00388 00389 // This rule is originally from Keast: 00390 // Patrick Keast, 00391 // Moderate Degree Tetrahedral Quadrature Formulas, 00392 // Computer Methods in Applied Mechanics and Engineering, 00393 // Volume 55, Number 3, May 1986, pages 339-348. 00394 // 00395 // It is accurate on 6th-degree polynomials and has 24 points 00396 // vs. 64 for the comparable conical product rule. 00397 // 00398 // Values copied 24th June 2008 from: 00399 // http://people.scs.fsu.edu/~burkardt/f_src/keast/keast.f90 00400 case SIXTH: 00401 { 00402 _points.resize (24); 00403 _weights.resize(24); 00404 00405 // The raw data for the quadrature rule. 00406 const Real rule_data[4][4] = { 00407 {0.356191386222544953e+00 , 0.214602871259151684e+00 , 0., 0.00665379170969464506e+00}, // 4 00408 {0.877978124396165982e+00 , 0.0406739585346113397e+00, 0., 0.00167953517588677620e+00}, // 4 00409 {0.0329863295731730594e+00, 0.322337890142275646e+00 , 0., 0.00922619692394239843e+00}, // 4 00410 {0.0636610018750175299e+00, 0.269672331458315867e+00 , 0.603005664791649076e+00, 0.00803571428571428248e+00} // 12 00411 }; 00412 00413 00414 // Now call the keast routine to generate _points and _weights 00415 keast_rule(rule_data, 4); 00416 00417 return; 00418 } 00419 00420 00421 // Keast's 31 point, 7th-order rule contains points on the reference 00422 // element boundary, so we've decided not to include it here. 00423 // 00424 // Keast's 8th-order rule has 45 points. and a negative 00425 // weight, so if you've explicitly disallowed such rules 00426 // you will fall through to the conical product rules 00427 // below. 00428 case SEVENTH: 00429 case EIGHTH: 00430 { 00431 if (allow_rules_with_negative_weights) 00432 { 00433 _points.resize (45); 00434 _weights.resize(45); 00435 00436 // The raw data for the quadrature rule. 00437 const Real rule_data[7][4] = { 00438 {0.250000000000000000e+00, 0., 0., -0.393270066412926145e-01}, // 1 00439 {0.617587190300082967e+00, 0.127470936566639015e+00, 0., 0.408131605934270525e-02}, // 4 00440 {0.903763508822103123e+00, 0.320788303926322960e-01, 0., 0.658086773304341943e-03}, // 4 00441 {0.497770956432810185e-01, 0., 0.450222904356718978e+00, 0.438425882512284693e-02}, // 6 00442 {0.183730447398549945e+00, 0., 0.316269552601450060e+00, 0.138300638425098166e-01}, // 6 00443 {0.231901089397150906e+00, 0.229177878448171174e-01, 0.513280033360881072e+00, 0.424043742468372453e-02}, // 12 00444 {0.379700484718286102e-01, 0.730313427807538396e+00, 0.193746475248804382e+00, 0.223873973961420164e-02} // 12 00445 }; 00446 00447 00448 // Now call the keast routine to generate _points and _weights 00449 keast_rule(rule_data, 7); 00450 00451 return; 00452 } // end if (allow_rules_with_negative_weights) 00453 // Note: if !allow_rules_with_negative_weights, fall through to next case. 00454 } 00455 00456 // Fall back on Grundmann-Moller or Conical Product rules at high orders. 00457 default: 00458 { 00459 if ((allow_rules_with_negative_weights) && (_order + 2*p < 34)) 00460 { 00461 // The Grundmann-Moller rules are defined to arbitrary order and 00462 // can have significantly fewer evaluation points than concial product 00463 // rules. If you allow rules with negative weights, the GM rules 00464 // will be more efficient for degree up to 33 (for degree 34 and 00465 // higher, CP is more efficient!) but may be more susceptible 00466 // to round-off error. Safest is to disallow rules with negative 00467 // weights, but this decision should be made on a case-by-case basis. 00468 QGrundmann_Moller gm_rule(3, _order); 00469 gm_rule.init(type_in, p); 00470 00471 // Swap points and weights with the about-to-be destroyed rule. 00472 _points.swap (gm_rule.get_points() ); 00473 _weights.swap(gm_rule.get_weights()); 00474 00475 return; 00476 } 00477 00478 else 00479 { 00480 // The following quadrature rules are generated as 00481 // conical products. These tend to be non-optimal 00482 // (use too many points, cluster points in certain 00483 // regions of the domain) but they are quite easy to 00484 // automatically generate using a 1D Gauss rule on 00485 // [0,1] and two 1D Jacobi-Gauss rules on [0,1]. 00486 QConical conical_rule(3, _order); 00487 conical_rule.init(type_in, p); 00488 00489 // Swap points and weights with the about-to-be destroyed rule. 00490 _points.swap (conical_rule.get_points() ); 00491 _weights.swap(conical_rule.get_weights()); 00492 00493 return; 00494 } 00495 } 00496 } 00497 } // end case TET4,TET10 00498 00499 00500 00501 //--------------------------------------------- 00502 // Prism quadrature rules 00503 case PRISM6: 00504 case PRISM15: 00505 case PRISM18: 00506 { 00507 // We compute the 3D quadrature rule as a tensor 00508 // product of the 1D quadrature rule and a 2D 00509 // triangle quadrature rule 00510 00511 QGauss q1D(1,_order); 00512 QGauss q2D(2,_order); 00513 00514 // Initialize 00515 q1D.init(EDGE2,p); 00516 q2D.init(TRI3,p); 00517 00518 tensor_product_prism(q1D, q2D); 00519 00520 return; 00521 } 00522 00523 00524 00525 //--------------------------------------------- 00526 // Pyramid 00527 case PYRAMID5: 00528 case PYRAMID13: 00529 case PYRAMID14: 00530 { 00531 // We compute the Pyramid rule as a conical product of a 00532 // Jacobi rule with alpha==2 on the interval [0,1] two 1D 00533 // Gauss rules with suitably modified points. The idea comes 00534 // from: Stroud, A.H. "Approximate Calculation of Multiple 00535 // Integrals." 00536 QConical conical_rule(3, _order); 00537 conical_rule.init(type_in, p); 00538 00539 // Swap points and weights with the about-to-be destroyed rule. 00540 _points.swap (conical_rule.get_points() ); 00541 _weights.swap(conical_rule.get_weights()); 00542 00543 return; 00544 00545 } 00546 00547 00548 00549 //--------------------------------------------- 00550 // Unsupported type 00551 default: 00552 libmesh_error_msg("ERROR: Unsupported type: " << type_in); 00553 } 00554 #endif 00555 } 00556 00557 } // namespace libMesh