Branch data Line data Source code
1 : : // ***************************************************************************** 2 : : /*! 3 : : \file src/PDE/Integrate/Quadrature.cpp 4 : : \copyright 2012-2015 J. Bakosi, 5 : : 2016-2018 Los Alamos National Security, LLC., 6 : : 2019-2021 Triad National Security, LLC. 7 : : All rights reserved. See the LICENSE file for details. 8 : : \brief Quadrature coordinates and weights for numerical integration 9 : : \details This file contains functions that provide Gauss quadrature 10 : : coordinates and weights for numerical integration. 11 : : */ 12 : : // ***************************************************************************** 13 : : 14 : : #include "Quadrature.hpp" 15 : : 16 : : void 17 : 26349595 : tk::GaussQuadratureTet( const std::size_t NG, 18 : : std::array< std::vector< real >, 3>& coordgp, 19 : : std::vector< real >& wgp ) 20 : : // ***************************************************************************** 21 : : //! Initialize Gaussian quadrature points locations and weights for a tetrahedron 22 : : //! \param[in] NG number of quadrature points 23 : : //! \param[in,out] coordgp 3 spatial coordinates of quadrature points 24 : : //! \param[in,out] wgp Weights of quadrature points 25 : : // ***************************************************************************** 26 : : { 27 [ - + ][ - - ]: 26349595 : Assert( coordgp[0].size() == NG, "Size mismatch" ); [ - - ][ - - ] 28 [ - + ][ - - ]: 26349595 : Assert( coordgp[1].size() == NG, "Size mismatch" ); [ - - ][ - - ] 29 [ - + ][ - - ]: 26349595 : Assert( coordgp[1].size() == NG, "Size mismatch" ); [ - - ][ - - ] 30 [ - + ][ - - ]: 26349595 : Assert( wgp.size() == NG, "Size mismatch" ); [ - - ][ - - ] 31 : : 32 [ + + ][ + + ]: 26349595 : switch( NG ) [ + - ] 33 : : { 34 : 15321251 : case 1: 35 : : { 36 : 15321251 : coordgp[0][0] = 0.25; 37 : 15321251 : coordgp[1][0] = 0.25; 38 : 15321251 : coordgp[2][0] = 0.25; 39 : 15321251 : wgp[0] = 1.0; 40 : : } 41 : 15321251 : break; 42 : : 43 : 457100 : case 4: 44 : : { 45 : 457100 : const real a1 = 0.5854101966249685; 46 : 457100 : const real a2 = 0.1381966011250105; 47 : : 48 : 457100 : coordgp[0][0] = a2; 49 : 457100 : coordgp[1][0] = a2; 50 : 457100 : coordgp[2][0] = a2; 51 : 457100 : wgp[0] = 0.25; 52 : : 53 : 457100 : coordgp[0][1] = a1; 54 : 457100 : coordgp[1][1] = a2; 55 : 457100 : coordgp[2][1] = a2; 56 : 457100 : wgp[1] = 0.25; 57 : : 58 : 457100 : coordgp[0][2] = a2; 59 : 457100 : coordgp[1][2] = a1; 60 : 457100 : coordgp[2][2] = a2; 61 : 457100 : wgp[2] = 0.25; 62 : : 63 : 457100 : coordgp[0][3] = a2; 64 : 457100 : coordgp[1][3] = a2; 65 : 457100 : coordgp[2][3] = a1; 66 : 457100 : wgp[3] = 0.25; 67 : : } 68 : 457100 : break; 69 : : 70 : 7048627 : case 5: 71 : : { 72 : 7048627 : coordgp[0][0] = 0.25; 73 : 7048627 : coordgp[1][0] = 0.25; 74 : 7048627 : coordgp[2][0] = 0.25; 75 : 7048627 : wgp[0] = -12.0/15.0; 76 : : 77 : 7048627 : coordgp[0][1] = 1.0/6.0; 78 : 7048627 : coordgp[1][1] = 1.0/6.0; 79 : 7048627 : coordgp[2][1] = 1.0/6.0; 80 : 7048627 : wgp[1] = 9.0/20.0; 81 : : 82 : 7048627 : coordgp[0][2] = 0.5; 83 : 7048627 : coordgp[1][2] = 1.0/6.0; 84 : 7048627 : coordgp[2][2] = 1.0/6.0; 85 : 7048627 : wgp[2] = 9.0/20.0; 86 : : 87 : 7048627 : coordgp[0][3] = 1.0/6.0; 88 : 7048627 : coordgp[1][3] = 0.5; 89 : 7048627 : coordgp[2][3] = 1.0/6.0; 90 : 7048627 : wgp[3] = 9.0/20.0; 91 : : 92 : 7048627 : coordgp[0][4] = 1.0/6.0; 93 : 7048627 : coordgp[1][4] = 1.0/6.0; 94 : 7048627 : coordgp[2][4] = 0.5; 95 : 7048627 : wgp[4] = 9.0/20.0; 96 : : } 97 : 7048627 : break; 98 : : 99 : 3411429 : case 11: 100 : : { 101 : 3411429 : const tk::real c1 = 0.3994035761667992; 102 : 3411429 : const tk::real c2 = 0.1005964238332008; 103 : 3411429 : const tk::real c3 = 343.0 / 7500.0; 104 : 3411429 : const tk::real c4 = 56.0 / 375.0; 105 : : 106 : 3411429 : coordgp[0][0] = 0.25; 107 : 3411429 : coordgp[1][0] = 0.25; 108 : 3411429 : coordgp[2][0] = 0.25; 109 : 3411429 : wgp[0] = -148.0/1875.0; 110 : : 111 : 3411429 : coordgp[0][1] = 11.0/14.0; 112 : 3411429 : coordgp[1][1] = 1.0/14.0; 113 : 3411429 : coordgp[2][1] = 1.0/14.0; 114 : 3411429 : wgp[1] = c3; 115 : : 116 : 3411429 : coordgp[0][2] = 1.0/14.0; 117 : 3411429 : coordgp[1][2] = 11.0/14.0; 118 : 3411429 : coordgp[2][2] = 1.0/14.0; 119 : 3411429 : wgp[2] = c3; 120 : : 121 : 3411429 : coordgp[0][3] = 1.0/14.0; 122 : 3411429 : coordgp[1][3] = 1.0/14.0; 123 : 3411429 : coordgp[2][3] = 11.0/14.0; 124 : 3411429 : wgp[3] = c3; 125 : : 126 : 3411429 : coordgp[0][4] = 1.0/14.0; 127 : 3411429 : coordgp[1][4] = 1.0/14.0; 128 : 3411429 : coordgp[2][4] = 1.0/14.0; 129 : 3411429 : wgp[4] = c3; 130 : : 131 : 3411429 : coordgp[0][5] = c1; 132 : 3411429 : coordgp[1][5] = c1; 133 : 3411429 : coordgp[2][5] = c2; 134 : 3411429 : wgp[5] = c4; 135 : : 136 : 3411429 : coordgp[0][6] = c1; 137 : 3411429 : coordgp[1][6] = c2; 138 : 3411429 : coordgp[2][6] = c1; 139 : 3411429 : wgp[6] = c4; 140 : : 141 : 3411429 : coordgp[0][7] = c1; 142 : 3411429 : coordgp[1][7] = c2; 143 : 3411429 : coordgp[2][7] = c2; 144 : 3411429 : wgp[7] = c4; 145 : : 146 : 3411429 : coordgp[0][8] = c2; 147 : 3411429 : coordgp[1][8] = c1; 148 : 3411429 : coordgp[2][8] = c1; 149 : 3411429 : wgp[8] = c4; 150 : : 151 : 3411429 : coordgp[0][9] = c2; 152 : 3411429 : coordgp[1][9] = c1; 153 : 3411429 : coordgp[2][9] = c2; 154 : 3411429 : wgp[9] = c4; 155 : : 156 : 3411429 : coordgp[0][10] = c2; 157 : 3411429 : coordgp[1][10] = c2; 158 : 3411429 : coordgp[2][10] = c1; 159 : 3411429 : wgp[10] = c4; 160 : : } 161 : 3411429 : break; 162 : : 163 : 111188 : case 14: 164 : : { 165 : 111188 : const tk::real a = 0.0673422422100983; 166 : 111188 : const tk::real b = 0.3108859192633005; 167 : 111188 : const tk::real c = 0.7217942490673264; 168 : 111188 : const tk::real d = 0.0927352503108912; 169 : 111188 : const tk::real e = 0.4544962958743506; 170 : 111188 : const tk::real f = 0.0455037041256494; 171 : 111188 : const tk::real p = 0.1126879257180162; 172 : 111188 : const tk::real q = 0.0734930431163619; 173 : 111188 : const tk::real r = 0.0425460207770812; 174 : : 175 : 111188 : coordgp[0][0] = a; 176 : 111188 : coordgp[1][0] = b; 177 : 111188 : coordgp[2][0] = b; 178 : 111188 : wgp[0] = p; 179 : : 180 : 111188 : coordgp[0][1] = b; 181 : 111188 : coordgp[1][1] = a; 182 : 111188 : coordgp[2][1] = b; 183 : 111188 : wgp[1] = p; 184 : : 185 : 111188 : coordgp[0][2] = b; 186 : 111188 : coordgp[1][2] = b; 187 : 111188 : coordgp[2][2] = a; 188 : 111188 : wgp[2] = p; 189 : : 190 : 111188 : coordgp[0][3] = b; 191 : 111188 : coordgp[1][3] = b; 192 : 111188 : coordgp[2][3] = b; 193 : 111188 : wgp[3] = p; 194 : : 195 : 111188 : coordgp[0][4] = c; 196 : 111188 : coordgp[1][4] = d; 197 : 111188 : coordgp[2][4] = d; 198 : 111188 : wgp[4] = q; 199 : : 200 : 111188 : coordgp[0][5] = d; 201 : 111188 : coordgp[1][5] = c; 202 : 111188 : coordgp[2][5] = d; 203 : 111188 : wgp[5] = q; 204 : : 205 : 111188 : coordgp[0][6] = d; 206 : 111188 : coordgp[1][6] = d; 207 : 111188 : coordgp[2][6] = c; 208 : 111188 : wgp[6] = q; 209 : : 210 : 111188 : coordgp[0][7] = d; 211 : 111188 : coordgp[1][7] = d; 212 : 111188 : coordgp[2][7] = d; 213 : 111188 : wgp[7] = q; 214 : : 215 : 111188 : coordgp[0][8] = e; 216 : 111188 : coordgp[1][8] = e; 217 : 111188 : coordgp[2][8] = f; 218 : 111188 : wgp[8] = r; 219 : : 220 : 111188 : coordgp[0][9] = e; 221 : 111188 : coordgp[1][9] = f; 222 : 111188 : coordgp[2][9] = e; 223 : 111188 : wgp[9] = r; 224 : : 225 : 111188 : coordgp[0][10] = e; 226 : 111188 : coordgp[1][10] = f; 227 : 111188 : coordgp[2][10] = f; 228 : 111188 : wgp[10] = r; 229 : : 230 : 111188 : coordgp[0][11] = f; 231 : 111188 : coordgp[1][11] = e; 232 : 111188 : coordgp[2][11] = e; 233 : 111188 : wgp[11] = r; 234 : : 235 : 111188 : coordgp[0][12] = f; 236 : 111188 : coordgp[1][12] = e; 237 : 111188 : coordgp[2][12] = f; 238 : 111188 : wgp[12] = r; 239 : : 240 : 111188 : coordgp[0][13] = f; 241 : 111188 : coordgp[1][13] = f; 242 : 111188 : coordgp[2][13] = e; 243 : 111188 : wgp[13] = r; 244 : : } 245 : 111188 : break; 246 : : } 247 : 26349595 : } 248 : : 249 : : void 250 : 88340128 : tk::GaussQuadratureTri( const std::size_t NG, 251 : : std::array< std::vector< real >, 2>& coordgp, 252 : : std::vector< real >& wgp ) 253 : : // ***************************************************************************** 254 : : //! Initialize Gaussian quadrature points locations and weights for a triangle 255 : : //! \param[in] NG number of quadrature points 256 : : //! \param[in,out] coordgp 2 spatial coordinates of quadrature points 257 : : //! \param[in,out] wgp Weights of quadrature points 258 : : // ***************************************************************************** 259 : : { 260 [ - + ][ - - ]: 88340128 : Assert( coordgp[0].size() == NG, "Size mismatch" ); [ - - ][ - - ] 261 [ - + ][ - - ]: 88340128 : Assert( coordgp[1].size() == NG, "Size mismatch" ); [ - - ][ - - ] 262 [ - + ][ - - ]: 88340128 : Assert( wgp.size() == NG, "Size mismatch" ); [ - - ][ - - ] 263 : : 264 [ + + ][ - + ]: 88340128 : switch( NG ) [ - ] 265 : : { 266 : 69268130 : case 1: 267 : 69268130 : coordgp[0][0] = 1.0/3.0; 268 : 69268130 : coordgp[1][0] = 1.0/3.0; 269 : 69268130 : wgp[0] = 1.0; 270 : 69268130 : break; 271 : : 272 : 14362228 : case 3: 273 : 14362228 : coordgp[0][0] = 2.0/3.0; 274 : 14362228 : coordgp[1][0] = 1.0/6.0; 275 : 14362228 : wgp[0] = 1.0/3.0; 276 : : 277 : 14362228 : coordgp[0][1] = 1.0/6.0; 278 : 14362228 : coordgp[1][1] = 2.0/3.0; 279 : 14362228 : wgp[1] = 1.0/3.0; 280 : : 281 : 14362228 : coordgp[0][2] = 1.0/6.0; 282 : 14362228 : coordgp[1][2] = 1.0/6.0; 283 : 14362228 : wgp[2] = 1.0/3.0; 284 : 14362228 : break; 285 : : 286 : 0 : case 4: 287 : 0 : coordgp[0][0] = 1.0/3.0; 288 : 0 : coordgp[1][0] = 1.0/3.0; 289 : 0 : wgp[0] = -27.0/48.0; 290 : : 291 : 0 : coordgp[0][1] = 1.0/5.0; 292 : 0 : coordgp[1][1] = 1.0/5.0; 293 : 0 : wgp[1] = 25.0/48.0; 294 : : 295 : 0 : coordgp[0][2] = 3.0/5.0; 296 : 0 : coordgp[1][2] = 1.0/5.0; 297 : 0 : wgp[2] = 25.0/48.0; 298 : : 299 : 0 : coordgp[0][3] = 1.0/5.0; 300 : 0 : coordgp[1][3] = 3.0/5.0; 301 : 0 : wgp[3] = 25.0/48.0; 302 : 0 : break; 303 : : 304 : 4709770 : case 6: 305 : 4709770 : const tk::real c1 = 0.816847572980459; 306 : 4709770 : const tk::real c2 = 0.091576213509771; 307 : 4709770 : const tk::real c3 = 0.091576213509771; 308 : 4709770 : const tk::real c4 = 0.108103018168070; 309 : 4709770 : const tk::real c5 = 0.445948490915965; 310 : 4709770 : const tk::real c6 = 0.445948490915965; 311 : 4709770 : const tk::real w1 = 0.054975870996713638 * 2.0; 312 : 4709770 : const tk::real w2 = 0.1116907969117165 * 2.0; 313 : : 314 : 4709770 : coordgp[0][0] = c1; 315 : 4709770 : coordgp[1][0] = c2; 316 : 4709770 : wgp[0] = w1; 317 : : 318 : 4709770 : coordgp[0][1] = c2; 319 : 4709770 : coordgp[1][1] = c3; 320 : 4709770 : wgp[1] = w1; 321 : : 322 : 4709770 : coordgp[0][2] = c3; 323 : 4709770 : coordgp[1][2] = c1; 324 : 4709770 : wgp[2] = w1; 325 : : 326 : 4709770 : coordgp[0][3] = c4; 327 : 4709770 : coordgp[1][3] = c5; 328 : 4709770 : wgp[3] = w2; 329 : : 330 : 4709770 : coordgp[0][4] = c5; 331 : 4709770 : coordgp[1][4] = c6; 332 : 4709770 : wgp[4] = w2; 333 : : 334 : 4709770 : coordgp[0][5] = c6; 335 : 4709770 : coordgp[1][5] = c4; 336 : 4709770 : wgp[5] = w2; 337 : 4709770 : break; 338 : : } 339 : 88340128 : }