Line data Source code
1 : // The libMesh Finite Element Library. 2 : // Copyright (C) 2002-2026 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_trap.h" 22 : #include "libmesh/enum_to_string.h" 23 : 24 : namespace libMesh 25 : { 26 : 27 1976 : void QTrap::init_3D() 28 : { 29 : #if LIBMESH_DIM == 3 30 : 31 : //----------------------------------------------------------------------- 32 : // 3D quadrature rules 33 1976 : switch (_type) 34 : { 35 : //--------------------------------------------- 36 : // Hex quadrature rules 37 782 : case HEX8: 38 : case HEX20: 39 : case HEX27: 40 : { 41 : // We compute the 3D quadrature rule as a tensor 42 : // product of the 1D quadrature rule. 43 978 : QTrap q1D(1); 44 : 45 782 : tensor_product_hex( q1D ); 46 : 47 196 : return; 48 : } 49 : 50 : 51 : 52 : //--------------------------------------------- 53 : // Tetrahedral quadrature rules 54 14 : case TET4: 55 : case TET10: 56 : case TET14: 57 : { 58 14 : _points.resize(4); 59 14 : _weights.resize(4); 60 : 61 14 : _points[0](0) = 0.; 62 14 : _points[0](1) = 0.; 63 14 : _points[0](2) = 0.; 64 : 65 14 : _points[1](0) = 1.; 66 14 : _points[1](1) = 0.; 67 14 : _points[1](2) = 0.; 68 : 69 14 : _points[2](0) = 0.; 70 14 : _points[2](1) = 1.; 71 14 : _points[2](2) = 0.; 72 : 73 14 : _points[3](0) = 0.; 74 14 : _points[3](1) = 0.; 75 14 : _points[3](2) = 1.; 76 : 77 : 78 : 79 14 : _weights[0] = 1/Real(24); 80 14 : _weights[1] = _weights[0]; 81 14 : _weights[2] = _weights[0]; 82 14 : _weights[3] = _weights[0]; 83 : 84 14 : return; 85 : } 86 : 87 : 88 : //--------------------------------------------- 89 : // Pyramid quadrature rules 90 782 : case PYRAMID5: 91 : case PYRAMID13: 92 : case PYRAMID14: 93 : case PYRAMID18: 94 : { 95 782 : libmesh_error_msg_if(!allow_nodal_pyramid_quadrature, 96 : "Nodal quadrature on Pyramid elements is not allowed by default since\n" 97 : "the Jacobian of the inverse element map is not well-defined at the Pyramid apex.\n" 98 : "Set the QBase::allow_nodal_pyramid_quadrature flag to true to ignore skip this check."); 99 : 100 782 : _points.resize(5); 101 782 : _weights.resize(5); 102 : 103 782 : _points[0](0) = -1.; 104 782 : _points[0](1) = -1.; 105 782 : _points[0](2) = 0.; 106 : 107 782 : _points[1](0) = 1.; 108 782 : _points[1](1) = -1.; 109 782 : _points[1](2) = 0.; 110 : 111 782 : _points[2](0) = 1.; 112 782 : _points[2](1) = 1.; 113 782 : _points[2](2) = 0.; 114 : 115 782 : _points[3](0) = -1.; 116 782 : _points[3](1) = 1.; 117 782 : _points[3](2) = 0.; 118 : 119 782 : _points[4](0) = 0.; 120 782 : _points[4](1) = 0.; 121 782 : _points[4](2) = 1.; 122 : 123 : 124 : // These are of dubious value since we can't integrate on the 125 : // vertex where the mapping Jacobian is ill-defined, but if we 126 : // could, this is what would give exact solutions for 127 : // constants and linears on the master element. 128 782 : _weights[0] = 1/Real(4); 129 782 : _weights[1] = _weights[0]; 130 782 : _weights[2] = _weights[0]; 131 782 : _weights[3] = _weights[0]; 132 782 : _weights[4] = 1/Real(3); 133 : 134 782 : return; 135 : } 136 : 137 : //--------------------------------------------- 138 : // Prism quadrature rules 139 398 : case PRISM6: 140 : case PRISM15: 141 : case PRISM18: 142 : case PRISM20: 143 : case PRISM21: 144 : { 145 : // We compute the 3D quadrature rule as a tensor 146 : // product of the 1D quadrature rule and a 2D 147 : // triangle quadrature rule 148 : 149 498 : QTrap q1D(1); 150 200 : QTrap q2D(2); 151 : 152 : // Initialize the 2D rule (1D is pre-initialized) 153 398 : q2D.init(TRI3, _p_level, /*simple_type_only=*/true); 154 : 155 398 : tensor_product_prism(q1D, q2D); 156 : 157 100 : return; 158 : } 159 : 160 : 161 : //--------------------------------------------- 162 : // Unsupported type 163 0 : default: 164 0 : libmesh_error_msg("ERROR: Unsupported type: " << Utility::enum_to_string(_type)); 165 : } 166 : #endif 167 : } 168 : 169 : } // namespace libMesh