20 #include "libmesh/quadrature_conical.h" 
   21 #include "libmesh/quadrature_gauss.h" 
   22 #include "libmesh/quadrature_jacobi.h" 
   23 #include "libmesh/enum_quadrature_type.h" 
   56   libmesh_assert_equal_to (this->
get_dim(), 2);
 
   62   std::pair<Real, Real> old_range(-1.0L, 1.0L);
 
   63   std::pair<Real, Real> new_range( 0.0L, 1.0L);
 
   64   gauss1D.
scale(old_range,
 
   73   const unsigned int np = gauss1D.
n_points();
 
   76   libmesh_assert_greater_equal (gauss1D.
qp(0)(0), 0.0);
 
   77   libmesh_assert_less_equal (gauss1D.
qp(np-1)(0), 1.0);
 
   78   libmesh_assert_greater_equal (jac1D.
qp(0)(0), 0.0);
 
   79   libmesh_assert_less_equal (jac1D.
qp(np-1)(0), 1.0);
 
   87   for (
unsigned int i=0; i<np; i++)
 
   88     for (
unsigned int j=0; j<np; j++)
 
   91         _points[gp](1) = gauss1D.
qp(i)(0) * (1.-jac1D.
qp(j)(0)); 
 
  107   libmesh_assert_equal_to (this->
get_dim(), 3);
 
  114   std::pair<Real, Real> old_range(-1.0L, 1.0L);
 
  115   std::pair<Real, Real> new_range( 0.0L, 1.0L);
 
  116   gauss1D.
scale(old_range,
 
  126   const unsigned int np = gauss1D.
n_points();
 
  129   libmesh_assert_greater_equal (gauss1D.
qp(0)(0), 0.0);
 
  130   libmesh_assert_less_equal (gauss1D.
qp(np-1)(0), 1.0);
 
  131   libmesh_assert_greater_equal (jacA1D.
qp(0)(0), 0.0);
 
  132   libmesh_assert_less_equal (jacA1D.
qp(np-1)(0), 1.0);
 
  133   libmesh_assert_greater_equal (jacB1D.
qp(0)(0), 0.0);
 
  134   libmesh_assert_less_equal (jacB1D.
qp(np-1)(0), 1.0);
 
  142   for (
unsigned int i=0; i<np; i++)
 
  143     for (
unsigned int j=0; j<np; j++)
 
  144       for (
unsigned int k=0; k<np; k++)
 
  147           _points[gp](1) = jacA1D.
qp(j)(0)  * (1.-jacB1D.
qp(k)(0));                         
 
  148           _points[gp](2) = gauss1D.
qp(i)(0) * (1.-jacA1D.
qp(j)(0)) * (1.-jacB1D.
qp(k)(0)); 
 
  149           _weights[gp]   = gauss1D.
w(i)     * jacA1D.
w(j)          * jacB1D.
w(k);          
 
  184   libmesh_assert_equal_to (this->
get_dim(), 3);
 
  193   const unsigned int np = gauss1D.
n_points();
 
  201   for (
unsigned int i=0; i<np; ++i)
 
  202     for (
unsigned int j=0; j<np; ++j)
 
  203       for (
unsigned int k=0; k<np; ++k, ++q)
 
  205           const Real xi=gauss1D.
qp(i)(0);
 
  206           const Real yj=gauss1D.
qp(j)(0);
 
  207           const Real zk=jac1D.
qp(k)(0);
 
  212           _weights[q]   = gauss1D.
w(i) * gauss1D.
w(j) * jac1D.
w(k);