21 #include "libmesh/quadrature_monomial.h"    22 #include "libmesh/quadrature_gauss.h"    54               const Real data[1][4] =
    59               const unsigned int rule_id[1] = {
    82               const Real sqrt19 = std::sqrt(
Real(19));            
    83               const Real tp     = std::sqrt(71440 + 6802*sqrt19); 
    88               const Real lambda =  std::sqrt(
Real(1919)/3285 - 148*sqrt19/3285 + 4*tp/3285);
    91               const Real xi     = -std::sqrt(
Real(1121)/3285 +  74*sqrt19/3285 - 2*tp/3285);
    94               const Real mu     =  std::sqrt(
Real(1121)/3285 +  74*sqrt19/3285 + 2*tp/3285);
    97               const Real gamma  =  std::sqrt(
Real(1919)/3285 - 148*sqrt19/3285 - 4*tp/3285);
   107               const Real B      = 
Real(1) / (
Real(260072)/133225  - 1520*sqrt19/133225 + (133 - 37*sqrt19)*tp/133225);
   110               const Real C      = 
Real(1) / (
Real(260072)/133225  - 1520*sqrt19/133225 - (133 - 37*sqrt19)*tp/133225);
   188                   const Real data[3][4] =
   190                       {0.00000000000000000000000000000000e+00_R, 0.00000000000000000000000000000000e+00_R, 0.00000000000000000000000000000000e+00_R, -1.27536231884057971014492753623188e+00_R},
   191                       {5.85540043769119907612630781744060e-01_R, 0.00000000000000000000000000000000e+00_R, 0.00000000000000000000000000000000e+00_R,  8.71111111111111111111111111111111e-01_R},
   192                       {6.94470135991704766602025803883310e-01_R, 9.37161638568208038511047377665396e-01_R, 4.15659267604065126239606672567031e-01_R,  1.68695652173913043478260869565217e-01_R}
   195                   const unsigned int rule_id[3] = {
   222                 r  = std::sqrt(
Real(6)/7),                                     
   223                 s  = std::sqrt((
Real(960) - 3*std::sqrt(
Real(28798))) / 2726), 
   224                 t  = std::sqrt((
Real(960) + 3*std::sqrt(
Real(28798))) / 2726), 
   225                 B1 = 
Real(8624)/29160,                                         
   226                 B2 = 
Real(2744)/29160,                                         
   227                 B3 = 8*(774*t*t - 230)/(9720*(t*t-s*s)),                       
   228                 B4 = 8*(230 - 774*s*s)/(9720*(t*t-s*s));                       
   230               const Real data[4][4] =
   238               const unsigned int rule_id[4] = {
   283               const Real data[5][4] =
   285                   {0.00000000000000000000000000000000e+00_R, 0.00000000000000000000000000000000e+00_R, 0.00000000000000000000000000000000e+00_R, 4.51903714875199690490763818699555e-01_R},
   286                   {7.82460796435951590652813975429717e-01_R, 0.00000000000000000000000000000000e+00_R, 0.00000000000000000000000000000000e+00_R, 2.99379177352338919703385618576171e-01_R},
   287                   {4.88094669706366480526729301468686e-01_R, 4.88094669706366480526729301468686e-01_R, 4.88094669706366480526729301468686e-01_R, 3.00876159371240019939698689791164e-01_R},
   288                   {8.62218927661481188856422891110042e-01_R, 8.62218927661481188856422891110042e-01_R, 8.62218927661481188856422891110042e-01_R, 4.94843255877038125738173175714853e-02_R},
   289                   {2.81113909408341856058098281846420e-01_R, 9.44196578292008195318687494773744e-01_R, 6.97574833707236996779391729948984e-01_R, 1.22872389222467338799199767122592e-01_R}
   292               const unsigned int rule_id[5] = {
   319       libmesh_fallthrough();
   325         gauss_rule.
init(*
this);
 
bool allow_rules_with_negative_weights
Flag (default true) controlling the use of quadrature rules with negative weights. 
ElemType _type
The type of element for which the current values have been computed. 
const std::vector< Real > & get_weights() const
The libMesh namespace provides an interface to certain functionality in the library. 
std::vector< Point > _points
The locations of the quadrature points in reference element space. 
std::vector< Real > _weights
The quadrature weights. 
void kim_rule(const Real rule_data[][4], const unsigned int *rule_id, const unsigned int n_pts)
Rules from Kim and Song, Comm. 
virtual void init_3D() override
Initializes the 3D quadrature rule by filling the points and weights vectors with the appropriate val...
Order _order
The polynomial order which the quadrature rule is capable of integrating exactly. ...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const std::vector< Point > & get_points() const
This class implements specific orders of Gauss quadrature. 
virtual void init(const Elem &e, unsigned int p_level=invalid_uint)
Initializes the data structures for a quadrature rule for the element e. 
A Point defines a location in LIBMESH_DIM dimensional Real space.