📄 quadrature_gauss_3d.c
字号:
// $Id: quadrature_gauss_3D.C 2901 2008-07-02 15:40:29Z jwpeterson $// The libMesh Finite Element Library.// Copyright (C) 2002-2007 Benjamin S. Kirk, John W. Peterson // This library is free software; you can redistribute it and/or// modify it under the terms of the GNU Lesser General Public// License as published by the Free Software Foundation; either// version 2.1 of the License, or (at your option) any later version. // This library is distributed in the hope that it will be useful,// but WITHOUT ANY WARRANTY; without even the implied warranty of// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU// Lesser General Public License for more details. // You should have received a copy of the GNU Lesser General Public// License along with this library; if not, write to the Free Software// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA// Local includes#include "quadrature_gauss.h"#include "quadrature_jacobi.h"#include "quadrature_gm.h"void QGauss::init_3D(const ElemType _type, unsigned int p){#if DIM == 3 //----------------------------------------------------------------------- // 3D quadrature rules switch (_type) { //--------------------------------------------- // Hex quadrature rules case HEX8: case HEX20: case HEX27: { // We compute the 3D quadrature rule as a tensor // product of the 1D quadrature rule. QGauss q1D(1,_order); q1D.init(EDGE2,p); tensor_product_hex( q1D ); return; } //--------------------------------------------- // Tetrahedral quadrature rules case TET4: case TET10: { // Taken from pg. 222 of "The finite element method," vol. 1 // ed. 5 by Zienkiewicz & Taylor switch(_order + 2*p) { case CONSTANT: case FIRST: { // Exact for linears _points.resize(1); _weights.resize(1); _points[0](0) = .25; _points[0](1) = .25; _points[0](2) = .25; _weights[0] = .1666666666666666666666666666666666666666666667; return; } case SECOND: { // Exact for quadratics _points.resize(4); _weights.resize(4); const Real a = .585410196624969; const Real b = .138196601125011; _points[0](0) = a; _points[0](1) = b; _points[0](2) = b; _points[1](0) = b; _points[1](1) = a; _points[1](2) = b; _points[2](0) = b; _points[2](1) = b; _points[2](2) = a; _points[3](0) = b; _points[3](1) = b; _points[3](2) = b; _weights[0] = .0416666666666666666666666666666666666666666667; _weights[1] = _weights[0]; _weights[2] = _weights[0]; _weights[3] = _weights[0]; return; } // Can be found in the class notes // http://www.cs.rpi.edu/~flaherje/FEM/fem6.ps // by Flaherty. // // Caution: this rule has a negative weight and may be // unsuitable for some problems. // Exact for cubics. // // Note: Keast (see ref. elsewhere in this file) also gives // a third-order rule with positive weights, but it contains points // on the ref. elt. boundary, making it less suitable for FEM calculations. case THIRD: { if (allow_rules_with_negative_weights) { _points.resize(5); _weights.resize(5); _points[0](0) = .25; _points[0](1) = .25; _points[0](2) = .25; _points[1](0) = .5; _points[1](1) = .16666666666666666666666666666666666666666667; _points[1](2) = .16666666666666666666666666666666666666666667; _points[2](0) = .16666666666666666666666666666666666666666667; _points[2](1) = .5; _points[2](2) = .16666666666666666666666666666666666666666667; _points[3](0) = .16666666666666666666666666666666666666666667; _points[3](1) = .16666666666666666666666666666666666666666667; _points[3](2) = .5; _points[4](0) = .16666666666666666666666666666666666666666667; _points[4](1) = .16666666666666666666666666666666666666666667; _points[4](2) = .16666666666666666666666666666666666666666667; _weights[0] = -.133333333333333333333333333333333333333333333; _weights[1] = .075; _weights[2] = _weights[1]; _weights[3] = _weights[1]; _weights[4] = _weights[1]; return; } // end if (allow_rules_with_negative_weights) // Note: if !allow_rules_with_negative_weights, fall through to next case. } // Originally a Keast rule, // Patrick Keast, // Moderate Degree Tetrahedral Quadrature Formulas, // Computer Methods in Applied Mechanics and Engineering, // Volume 55, Number 3, May 1986, pages 339-348. // // Can also be found the class notes // http://www.cs.rpi.edu/~flaherje/FEM/fem6.ps // by Flaherty. // // Caution: this rule has a negative weight and may be // unsuitable for some problems. case FOURTH: { if (allow_rules_with_negative_weights) { _points.resize(11); _weights.resize(11); // The raw data for the quadrature rule. const Real p[3][4] = { {0.250000000000000000e+00, 0., 0., -0.131555555555555556e-01}, // 1 {0.785714285714285714e+00, 0.714285714285714285e-01, 0., 0.762222222222222222e-02}, // 4 {0.399403576166799219e+00, 0., 0.100596423833200785e+00, 0.248888888888888889e-01} // 6 }; // Now call the keast routine to generate _points and _weights keast_rule(p, 3); return; } // end if (allow_rules_with_negative_weights) // Note: if !allow_rules_with_negative_weights, fall through to next case. } // Walkington's fifth-order 14-point rule from // "Quadrature on Simplices of Arbitrary Dimension" // // We originally had a Keast rule here, but this rule had // more points than an equivalent rule by Walkington and // also contained points on the boundary of the ref. elt, // making it less suitable for FEM calculations. case FIFTH: { _points.resize(14); _weights.resize(14); // permutations of these points and suitably-modified versions of // these points are the quadrature point locations const Real a[3] = {0.31088591926330060980, // a1 from the paper 0.092735250310891226402, // a2 from the paper 0.045503704125649649492}; // a3 from the paper // weights. a[] and w[] are the only floating-point inputs required // for this rule. const Real w[3] = {0.018781320953002641800, // w1 from the paper 0.012248840519393658257, // w2 from the paper 0.0070910034628469110730}; // w3 from the paper // The first two sets of 4 points are formed in a similar manner for (unsigned int i=0; i<2; ++i) { // Where we will insert values into _points and _weights const unsigned int offset=4*i; // Stuff points and weights values into their arrays const Real b = 1. - 3.*a[i]; // Here are the permutations. Order of these is not important, // all have the same weight _points[offset + 0] = Point(a[i], a[i], a[i]); _points[offset + 1] = Point(a[i], b, a[i]); _points[offset + 2] = Point( b, a[i], a[i]); _points[offset + 3] = Point(a[i], a[i], b); // These 4 points all have the same weights for (unsigned int j=0; j<4; ++j) _weights[offset + j] = w[i]; } // end for { // The third set contains 6 points and is formed a little differently const unsigned int offset = 8; const Real b = 0.5*(1. - 2.*a[2]); // Here are the permutations. Order of these is not important, // all have the same weight _points[offset + 0] = Point(b , b, a[2]); _points[offset + 1] = Point(b , a[2], a[2]); _points[offset + 2] = Point(a[2], a[2], b); _points[offset + 3] = Point(a[2], b, a[2]); _points[offset + 4] = Point( b, a[2], b); _points[offset + 5] = Point(a[2], b, b); // These 6 points all have the same weights for (unsigned int j=0; j<6; ++j) _weights[offset + j] = w[2]; } // Original rule by Keast, unsuitable because it has points on the // reference element boundary! // _points.resize(15); // _weights.resize(15); // _points[0](0) = 0.25; // _points[0](1) = 0.25; // _points[0](2) = 0.25; // { // const Real a = 0.; // const Real b = 0.333333333333333333333333333333333333333; // _points[1](0) = a; // _points[1](1) = b; // _points[1](2) = b; // _points[2](0) = b; // _points[2](1) = a; // _points[2](2) = b; // _points[3](0) = b; // _points[3](1) = b; // _points[3](2) = a; // _points[4](0) = b; // _points[4](1) = b; // _points[4](2) = b; // } // { // const Real a = 0.7272727272727272727272727272727272727272727272727272727;
⌨️ 快捷键说明
复制代码
Ctrl + C
搜索代码
Ctrl + F
全屏模式
F11
切换主题
Ctrl + Shift + D
显示快捷键
?
增大字号
Ctrl + =
减小字号
Ctrl + -