Line data Source code
1 : // The libMesh Finite Element Library. 2 : // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 3 : 4 : // This library is free software; you can redistribute it and/or 5 : // modify it under the terms of the GNU Lesser General Public 6 : // License as published by the Free Software Foundation; either 7 : // version 2.1 of the License, or (at your option) any later version. 8 : 9 : // This library is distributed in the hope that it will be useful, 10 : // but WITHOUT ANY WARRANTY; without even the implied warranty of 11 : // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 12 : // Lesser General Public License for more details. 13 : 14 : // You should have received a copy of the GNU Lesser General Public 15 : // License along with this library; if not, write to the Free Software 16 : // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 17 : 18 : 19 : 20 : // Local includes 21 : #include "libmesh/quadrature_grid.h" 22 : #include "libmesh/enum_to_string.h" 23 : 24 : namespace libMesh 25 : { 26 : 27 1278 : void QGrid::init_2D() 28 : { 29 : #if LIBMESH_DIM > 1 30 : 31 1278 : switch (_type) 32 : { 33 : 34 : 35 : //--------------------------------------------- 36 : // Quadrilateral quadrature rules 37 639 : case QUAD4: 38 : case QUADSHELL4: 39 : case QUAD8: 40 : case QUADSHELL8: 41 : case QUAD9: 42 : case QUADSHELL9: 43 : { 44 : // We compute the 2D quadrature rule as a tensor 45 : // product of the 1D quadrature rule. 46 : // 47 : // We ignore p_level in QGrid, since the "order" here is a 48 : // user-requested number of grid points, nothing to do with 49 : // the order of exactly-integrated polynomial spaces 50 657 : QGrid q1D(1, _order); 51 639 : tensor_product_quad( q1D ); 52 18 : return; 53 : } 54 : 55 : 56 : //--------------------------------------------- 57 : // Triangle quadrature rules 58 639 : case TRI3: 59 : case TRISHELL3: 60 : case TRI3SUBDIVISION: 61 : case TRI6: 62 : case TRI7: 63 : { 64 639 : const unsigned int np = (_order + 1)*(_order + 2)/2; 65 639 : const Real weight = Real(0.5)/np; 66 639 : const Real dx = Real(1)/(_order+1); 67 639 : _points.resize(np); 68 639 : _weights.resize(np); 69 : 70 18 : unsigned int pt = 0; 71 4491 : for (int i = 0; i != _order + 1; ++i) 72 : { 73 19383 : for (int j = 0; j != _order + 1 - i; ++j) 74 : { 75 15549 : _points[pt](0) = (i+0.5)*dx; 76 15549 : _points[pt](1) = (j+0.5)*dx; 77 15549 : _weights[pt] = weight; 78 15549 : pt++; 79 : } 80 : } 81 18 : return; 82 : } 83 : 84 : //--------------------------------------------- 85 : // Unsupported type 86 0 : default: 87 0 : libmesh_error_msg("Element type not supported!:" << Utility::enum_to_string(_type)); 88 : } 89 : #endif 90 : } 91 : 92 : } // namespace libMesh