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 : 28 1917 : void QGrid::init_3D() 29 : { 30 : #if LIBMESH_DIM == 3 31 : 32 1917 : switch (_type) 33 : { 34 : //--------------------------------------------- 35 : // Hex quadrature rules 36 639 : case HEX8: 37 : case HEX20: 38 : case HEX27: 39 : { 40 : // We compute the 3D quadrature rule as a tensor 41 : // product of the 1D quadrature rule. 42 : // 43 : // We ignore p_level in QGrid, since the "order" here is a 44 : // user-requested number of grid points, nothing to do with 45 : // the order of exactly-integrated polynomial spaces 46 657 : QGrid q1D(1,_order); 47 : 48 639 : tensor_product_hex( q1D ); 49 : 50 18 : return; 51 : } 52 : 53 : 54 : 55 : //--------------------------------------------- 56 : // Tetrahedral quadrature rules 57 639 : case TET4: 58 : case TET10: 59 : case TET14: 60 : { 61 657 : const unsigned int np = (_order+1)*(_order+2)*(_order+3)/6; 62 : // Master tet has 1x1 triangle base, height 1, so volume = 1/6 63 639 : const Real weight = Real(1)/Real(6)/np; 64 639 : const Real dx = Real(1)/(_order+1); 65 639 : _points.resize(np); 66 639 : _weights.resize(np); 67 : 68 18 : unsigned int pt = 0; 69 4491 : for (int i = 0; i != _order + 1; ++i) 70 : { 71 19383 : for (int j = 0; j != _order + 1 - i; ++j) 72 : { 73 66243 : for (int k = 0; k != _order + 1 - i - j; ++k) 74 : { 75 50694 : _points[pt](0) = (i+0.5)*dx; 76 50694 : _points[pt](1) = (j+0.5)*dx; 77 50694 : _points[pt](2) = (k+0.5)*dx; 78 50694 : _weights[pt] = weight; 79 50694 : pt++; 80 : } 81 : } 82 : } 83 18 : return; 84 : } 85 : 86 : 87 : // Prism quadrature rules 88 0 : case PRISM6: 89 : case PRISM15: 90 : case PRISM18: 91 : case PRISM20: 92 : case PRISM21: 93 : { 94 : // We compute the 3D quadrature rule as a tensor 95 : // product of the 1D quadrature rule and a 2D 96 : // triangle quadrature rule 97 : // 98 : // We ignore p_level in QGrid, since the "order" here is a 99 : // user-requested number of grid points, nothing to do with 100 : // the order of exactly-integrated polynomial spaces 101 0 : QGrid q1D(1,_order); 102 0 : QGrid q2D(2,_order); 103 : 104 : // Initialize the 2D rule (1D is pre-initialized) 105 0 : q2D.init(TRI3, 0, /*simple_type_only=*/true); 106 : 107 0 : tensor_product_prism(q1D, q2D); 108 : 109 0 : return; 110 : } 111 : 112 : 113 : 114 : //--------------------------------------------- 115 : // Pyramid 116 639 : case PYRAMID5: 117 : case PYRAMID13: 118 : case PYRAMID14: 119 : case PYRAMID18: 120 : { 121 657 : const unsigned int np = (_order+1)*(_order+2)*(_order+3)/6; 122 639 : _points.resize(np); 123 639 : _weights.resize(np); 124 : // Master pyramid has 2x2 base, height 1, so volume = 4/3 125 639 : const Real weight = Real(4)/Real(3)/np; 126 639 : const Real dx = Real(2)/(_order+1); 127 639 : const Real dz = Real(1)/(_order+1); 128 : 129 18 : unsigned int pt = 0; 130 4473 : for (int k = 0; k != _order + 1; ++k) 131 : { 132 19383 : for (int i = 0; i != _order + 1 - k; ++i) 133 : { 134 66243 : for (int j = 0; j != _order + 1 - k - i; ++j) 135 : { 136 53550 : _points[pt](0) = (i+0.5)*dx-1.0 + 137 50694 : (k+0.5)*dz; 138 50694 : _points[pt](1) = (j+0.5)*dx-1.0 + 139 1428 : (k+0.5)*dz; 140 50694 : _points[pt](2) = (k+0.5)*dz; 141 50694 : _weights[pt] = weight; 142 50694 : pt++; 143 : } 144 : } 145 : } 146 18 : return; 147 : } 148 : 149 : 150 : 151 : //--------------------------------------------- 152 : // Unsupported type 153 0 : default: 154 0 : libmesh_error_msg("ERROR: Unsupported type: " << Utility::enum_to_string(_type)); 155 : } 156 : #endif 157 : } 158 : 159 : } // namespace libMesh