dune-grid  2.1.1
quadraturerules.hh
Go to the documentation of this file.
00001 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
00002 // vi: set et ts=4 sw=2 sts=2:
00003 
00004 #ifndef DUNE_QUADRATURERULES_HH
00005 #define DUNE_QUADRATURERULES_HH
00006 
00007 #include<iostream>
00008 #include <limits>
00009 #include<vector>
00010 #include<map>
00011 
00012 #include<dune/common/fvector.hh>
00013 #include<dune/common/exceptions.hh>
00014 #include<dune/common/stdstreams.hh>
00015 #include<dune/common/geometrytype.hh>
00016 
00022 namespace Dune {
00023 
00028   class QuadratureOrderOutOfRange : public NotImplemented {};
00029   
00033   template<typename ct, int dim>
00034   class QuadraturePoint {
00035   public:
00036         // compile time parameters
00037         enum { d=dim };
00038         typedef ct CoordType;
00039 
00040     static const unsigned int dimension = d;
00041     typedef ct Field;
00042     typedef Dune::FieldVector<ct,dim> Vector;
00043 
00045         QuadraturePoint (const Vector& x, ct w) : local(x)
00046       {
00047         wght = w;
00048       }
00049 
00051         const Vector& position () const
00052       {
00053         return local;
00054       }
00055 
00057         const ct &weight () const
00058       {
00059         return wght;
00060       }
00061     
00062   protected:
00063         FieldVector<ct, dim> local;
00064         ct wght;
00065   };
00066 
00070   namespace QuadratureType {
00071     enum Enum {
00072       Gauss      = 0,
00073       
00074       Jacobian_1_0 = 1,
00075       Jacobian_2_0 = 2,
00076       
00077       Simpson    = 3,
00078       Trap       = 4,
00079       Grid       = 5,
00080       
00081       Clough     = 21,
00082       
00083       Invalid_Rule = 127
00084     };
00085   }
00086 
00090   template<typename ct, int dim>
00091   class QuadratureRule : public std::vector<QuadraturePoint<ct,dim> >
00092   {
00093   public:
00094 
00096     QuadratureRule() : delivered_order(-1) {}
00097 
00099     QuadratureRule(GeometryType t) : geometry_type(t), delivered_order(-1) {}
00100 
00102     QuadratureRule(GeometryType t, int order) : geometry_type(t), delivered_order(order) {}
00103 
00105         enum { d=dim };
00106 
00108         typedef ct CoordType;
00109 
00111         virtual int order () const { return delivered_order; }
00112 
00114         virtual GeometryType type () const { return geometry_type; }
00115     virtual ~QuadratureRule(){}
00116 
00119     typedef typename std::vector<QuadraturePoint<ct,dim> >::const_iterator iterator;
00120 
00121   protected:
00122     GeometryType geometry_type;
00123         int delivered_order;
00124     
00125     void tensor_product (const QuadratureRule<ct,1> & q1d)
00126       {
00127         // fill in all the gauss points
00128         int m = q1d.size();
00129         int n = power(m,dim);
00130         for (int i=0; i<n; i++)
00131                 {
00132                   // compute multi index for Gauss point
00133                   int x[dim];
00134                   int z = i;
00135                   for (int k=0; k<dim; ++k)
00136           {
00137             x[k] = z%m;
00138             z = z/m;
00139           }
00140 
00141                   // compute coordinates and weight
00142                   double weight = 1.0;
00143                   FieldVector<ct, dim> local;
00144                   for (int k=0; k<dim; k++) 
00145           {
00146             local[k] = q1d[x[k]].position()[0];
00147             weight *= q1d[x[k]].weight();
00148           }
00149 
00150                   // put in container
00151                   this->push_back(QuadraturePoint<ct,dim>(local,weight));
00152                 }
00153       }
00154 
00155         int power (int y, int d)
00156       {
00157         int m=1;
00158         for (int i=0; i<d; i++) m *= y;
00159         return m;
00160       }
00161   };
00162 
00163   // Forward declaration of the factory class,
00164   // needed internally by the QuadratureRules container class.
00165   template<typename ctype, int dim> class QuadratureRuleFactory;
00166 
00170   template<typename ctype, int dim>
00171   class QuadratureRules {
00172 
00174     typedef std::pair<GeometryType,int> QuadratureRuleKey;
00175 
00177     typedef Dune::QuadratureRule<ctype, dim> QuadratureRule;
00178 
00180     const QuadratureRule& _rule(const GeometryType& t, int p, QuadratureType::Enum qt=QuadratureType::Gauss)
00181       {
00182         static std::map<QuadratureRuleKey, QuadratureRule> _quadratureMap;
00183         QuadratureRuleKey key(t,p);
00184         if (_quadratureMap.find(key) == _quadratureMap.end()) {
00185           /*
00186             The rule must be acquired before we can store it.
00187             If we write this in one command, an invalid rule
00188             would get stored in case of an exception.
00189           */
00190           QuadratureRule rule =
00191             QuadratureRuleFactory<ctype,dim>::rule(t,p,qt);
00192           _quadratureMap[key] = rule;
00193         }
00194         return _quadratureMap[key];
00195       }
00197     static QuadratureRules& instance()
00198       {
00199         static QuadratureRules instance;
00200         return instance;
00201       }
00203     QuadratureRules () {};
00204   public:
00206     static const QuadratureRule& rule(const GeometryType& t, int p, QuadratureType::Enum qt=QuadratureType::Gauss)
00207       {
00208         return instance()._rule(t,p,qt);
00209       }
00211     static const QuadratureRule& rule(const GeometryType::BasicType t, int p, QuadratureType::Enum qt=QuadratureType::Gauss)
00212       {
00213         GeometryType gt(t,dim);
00214         return instance()._rule(gt,p,qt);
00215       }
00216   };
00217 
00218   /***********************************************************/
00219 
00223   template<typename ct, int dim>
00224   class CubeQuadratureRule :
00225     public QuadratureRule<ct,dim>
00226   {
00227   public:
00229         enum { d=dim };
00230 
00232         enum { highest_order=CubeQuadratureRule<ct,1>::highest_order };
00233 
00235         typedef ct CoordType;
00236 
00238         typedef CubeQuadratureRule value_type;
00239 
00240     ~CubeQuadratureRule(){}
00241   private:
00242     friend class QuadratureRuleFactory<ct,dim>;
00244         CubeQuadratureRule (int p) : QuadratureRule<ct,dim>(GeometryType(GeometryType::cube, d))
00245       {
00246         QuadratureRule<ct,1> q1D = QuadratureRules<ct,1>::rule(GeometryType::cube, p);
00247         this->tensor_product( q1D );
00248         this->delivered_order = q1D.order();
00249       }
00250 
00251   };
00252 
00255   template<typename ct>
00256   class CubeQuadratureRule<ct,0> :
00257     public QuadratureRule<ct,0>
00258   {
00259   public:
00260         // compile time parameters
00261         enum { d=0 };
00262         enum { dim=0 };
00263         enum { highest_order=61 };
00264         typedef ct CoordType;
00265         typedef CubeQuadratureRule value_type;
00266 
00267     ~CubeQuadratureRule(){}
00268   private:
00269     friend class QuadratureRuleFactory<ct,dim>;
00270     CubeQuadratureRule (int p):
00271       QuadratureRule<ct,0>(GeometryType(GeometryType::cube, 0))
00272     {
00273       FieldVector<ct, dim> point(0.0);
00274 
00275       if (p > highest_order)
00276         DUNE_THROW(QuadratureOrderOutOfRange, "Quadrature rule " << p << " not supported!");
00277 
00278       this->delivered_order = highest_order;
00279       this->push_back(QuadraturePoint<ct,dim>(point, 1.0));
00280     }
00281   };
00282  
00284   template<typename ct,
00285            bool fundamental = std::numeric_limits<ct>::is_specialized>
00286   struct CubeQuadratureInitHelper;
00287   template<typename ct>
00288   struct CubeQuadratureInitHelper<ct, true> {
00289     static void init(int p,
00290                      std::vector< FieldVector<ct, 1> > & _points,
00291                      std::vector< ct > & _weight,
00292                      int & delivered_order);
00293   };
00294   template<typename ct>
00295   struct CubeQuadratureInitHelper<ct, false> {
00296     static void init(int p,
00297                      std::vector< FieldVector<ct, 1> > & _points,
00298                      std::vector< ct > & _weight,
00299                      int & delivered_order);
00300   };
00301 
00304   template<typename ct>
00305   class CubeQuadratureRule<ct,1> :
00306     public QuadratureRule<ct,1>
00307   {
00308   public:
00309         // compile time parameters
00310         enum { d=1 };
00311         enum { dim=1 };
00312         enum { highest_order=61 };
00313         typedef ct CoordType;
00314         typedef CubeQuadratureRule value_type;
00315 
00316     ~CubeQuadratureRule(){}
00317   private:
00318     friend class QuadratureRuleFactory<ct,dim>;
00319     CubeQuadratureRule (int p)
00320       : QuadratureRule<ct,1>(GeometryType(GeometryType::cube, 1))
00321     {
00323       std::vector< FieldVector<ct, dim> > _points;
00324       std::vector< ct > _weight;
00325       
00326       CubeQuadratureInitHelper<ct>::init
00327         (p, _points, _weight, this->delivered_order);
00328       
00329       assert(_points.size() == _weight.size());
00330       for (size_t i = 0; i < _points.size(); i++)
00331         this->push_back(QuadraturePoint<ct,dim>(_points[i], _weight[i]));
00332     }
00333   };
00334 
00335   extern template CubeQuadratureRule<float, 1>::CubeQuadratureRule(int);
00336   extern template CubeQuadratureRule<double, 1>::CubeQuadratureRule(int);
00337 
00338 } // namespace Dune
00339 
00340 #define DUNE_INCLUDING_IMPLEMENTATION
00341 #include "quadraturerules/cube_imp.hh"
00342 
00343 namespace Dune {
00344 
00348   template<typename ct, int dim>
00349   class Jacobi1QuadratureRule;
00350   
00352   template<typename ct,
00353            bool fundamental = std::numeric_limits<ct>::is_specialized>
00354   struct Jacobi1QuadratureInitHelper;
00355   template<typename ct>
00356   struct Jacobi1QuadratureInitHelper<ct, true> {
00357     static void init(int p,
00358                      std::vector< FieldVector<ct, 1> > & _points,
00359                      std::vector< ct > & _weight,
00360                      int & delivered_order);
00361   };
00362   template<typename ct>
00363   struct Jacobi1QuadratureInitHelper<ct, false> {
00364     static void init(int p,
00365                      std::vector< FieldVector<ct, 1> > & _points,
00366                      std::vector< ct > & _weight,
00367                      int & delivered_order);
00368   };
00369 
00373   template<typename ct>
00374   class Jacobi1QuadratureRule<ct,1> :
00375     public QuadratureRule<ct,1>
00376   {
00377   public:
00379         enum { d=1 };
00382         enum { dim=1 };
00383 
00385         enum { highest_order=61 };
00386 
00388         typedef ct CoordType;
00389 
00391         typedef Jacobi1QuadratureRule value_type;
00392 
00393     ~Jacobi1QuadratureRule(){}
00394   private:
00395     friend class QuadratureRuleFactory<ct,dim>;
00396     Jacobi1QuadratureRule (int p)
00397       : QuadratureRule<ct,1>(GeometryType(GeometryType::cube, 1))
00398     {
00400       std::vector< FieldVector<ct, dim> > _points;
00401       std::vector< ct > _weight;
00402       
00403       int delivered_order;
00404       
00405       Jacobi1QuadratureInitHelper<ct>::init
00406         (p, _points, _weight, delivered_order);
00407       this->delivered_order = delivered_order;
00408       assert(_points.size() == _weight.size());
00409       for (size_t i = 0; i < _points.size(); i++)
00410         this->push_back(QuadraturePoint<ct,dim>(_points[i], _weight[i]));
00411     }
00412   };
00413 
00414 #ifndef DOXYGEN
00415   extern template Jacobi1QuadratureRule<float, 1>::Jacobi1QuadratureRule(int);
00416   extern template Jacobi1QuadratureRule<double, 1>::Jacobi1QuadratureRule(int);
00417 #endif // !DOXYGEN
00418 
00419 } // namespace Dune
00420 
00421 #define DUNE_INCLUDING_IMPLEMENTATION
00422 #include "quadraturerules/jacobi_1_0_imp.hh"
00423 
00424 namespace Dune {
00425  
00429   template<typename ct, int dim>
00430   class Jacobi2QuadratureRule;
00431   
00433   template<typename ct,
00434            bool fundamental = std::numeric_limits<ct>::is_specialized>
00435   struct Jacobi2QuadratureInitHelper;
00436   template<typename ct>
00437   struct Jacobi2QuadratureInitHelper<ct, true> {
00438     static void init(int p,
00439                      std::vector< FieldVector<ct, 1> > & _points,
00440                      std::vector< ct > & _weight,
00441                      int & delivered_order);
00442   };
00443   template<typename ct>
00444   struct Jacobi2QuadratureInitHelper<ct, false> {
00445     static void init(int p,
00446                      std::vector< FieldVector<ct, 1> > & _points,
00447                      std::vector< ct > & _weight,
00448                      int & delivered_order);
00449   };
00450 
00454   template<typename ct>
00455   class Jacobi2QuadratureRule<ct,1> :
00456     public QuadratureRule<ct,1>
00457   {
00458   public:
00460         enum { d=1 };
00461 
00464         enum { dim=1 };
00465 
00467         enum { highest_order=61 };
00468 
00470         typedef ct CoordType;
00471 
00473         typedef Jacobi2QuadratureRule value_type;
00474 
00475     ~Jacobi2QuadratureRule(){}
00476   private:
00477     friend class QuadratureRuleFactory<ct,dim>;
00478     Jacobi2QuadratureRule (int p)
00479       : QuadratureRule<ct,1>(GeometryType(GeometryType::cube, 1))
00480     {
00482       std::vector< FieldVector<ct, dim> > _points;
00483       std::vector< ct > _weight;
00484       
00485       int delivered_order;
00486       
00487       Jacobi2QuadratureInitHelper<ct>::init
00488         (p, _points, _weight, delivered_order);
00489       
00490       this->delivered_order = delivered_order;
00491       assert(_points.size() == _weight.size());
00492       for (size_t i = 0; i < _points.size(); i++)
00493         this->push_back(QuadraturePoint<ct,dim>(_points[i], _weight[i]));
00494     }
00495   };
00496 
00497 #ifndef DOXYGEN
00498   extern template Jacobi2QuadratureRule<float, 1>::Jacobi2QuadratureRule(int);
00499   extern template Jacobi2QuadratureRule<double, 1>::Jacobi2QuadratureRule(int);
00500 #endif // !DOXYGEN
00501 
00502 } // namespace Dune
00503 
00504 #define DUNE_INCLUDING_IMPLEMENTATION
00505 #include "quadraturerules/jacobi_2_0_imp.hh"
00506 
00507 namespace Dune {
00508  
00509   /************************************************
00510    * Quadraturerule for Simplices/Triangle
00511    *************************************************/
00512 
00516   template<typename ct, int dim>
00517   class SimplexQuadratureRule;
00518 
00522   template<typename ct>
00523   class SimplexQuadratureRule<ct,2> : public QuadratureRule<ct,2>
00524   {
00525   public:
00526 
00528     enum{d=2};
00529 
00531     enum { highest_order=CubeQuadratureRule<ct,1>::highest_order -1 };
00532 
00534     typedef ct CoordType;
00535 
00537     typedef SimplexQuadratureRule value_type;
00538     ~SimplexQuadratureRule(){}
00539   private:
00540     friend class QuadratureRuleFactory<ct,d>;
00541     SimplexQuadratureRule (int p);
00542   };
00543 
00547   template<typename ct>
00548   class SimplexQuadratureRule<ct,3> : public QuadratureRule<ct,3>
00549   {
00550   public:
00551 
00553     enum{d=3};
00554 
00556     enum { highest_order=CubeQuadratureRule<ct,1>::highest_order -2 };
00557 
00559     typedef ct CoordType;
00560 
00562     typedef SimplexQuadratureRule<ct,3> value_type;
00563     ~SimplexQuadratureRule(){}
00564   private:
00565     friend class QuadratureRuleFactory<ct,d>;
00566     SimplexQuadratureRule (int p);
00567   };
00568 
00569 /***********************************
00570  * quadrature for Prism
00571  **********************************/
00572 
00574   template<int dim>
00575   class PrismQuadraturePoints;
00576 
00578   template<>
00579   class PrismQuadraturePoints<3>
00580   {
00581   public:
00582         enum { MAXP=6};
00583         enum { highest_order=2 };
00584 
00586         PrismQuadraturePoints ()
00587       {
00588         int m = 0;
00589         O[m] = 0;
00590           
00591         // polynom degree 0  ???
00592         m = 6;
00593         G[m][0][0] = 0.0;
00594         G[m][0][1] = 0.0;
00595         G[m][0][2] = 0.0;
00596 
00597         G[m][1][0] = 1.0;
00598         G[m][1][1] = 0.0;
00599         G[m][1][2] = 0.0;
00600 
00601         G[m][2][0] = 0.0;
00602         G[m][2][1] = 1.0;
00603         G[m][2][2] = 0.0;
00604 
00605         G[m][3][0] = 0.0;
00606         G[m][3][1] = 0.0;
00607         G[m][3][2] = 1.0;
00608           
00609         G[m][4][0] = 1.0;
00610         G[m][4][1] = 0.0;
00611         G[m][4][2] = 1.0;
00612 
00613         G[m][5][0] = 0.0;
00614         G[m][5][1] = 0.1;
00615         G[m][5][2] = 1.0;
00616 
00617         W[m][0] = 0.16666666666666666 / 2.0;
00618         W[m][1] = 0.16666666666666666 / 2.0;
00619         W[m][2] = 0.16666666666666666 / 2.0;
00620         W[m][3] = 0.16666666666666666 / 2.0;
00621         W[m][4] = 0.16666666666666666 / 2.0;
00622         W[m][5] = 0.16666666666666666 / 2.0;
00623           
00624         O[m] = 0;// verify ????????
00625           
00626 
00627         // polynom degree 2  ???
00628         m = 6;
00629         G[m][0][0] =0.66666666666666666 ;
00630         G[m][0][1] =0.16666666666666666 ;
00631         G[m][0][2] =0.211324865405187 ;
00632 
00633         G[m][1][0] = 0.16666666666666666;
00634         G[m][1][1] =0.66666666666666666 ;
00635         G[m][1][2] = 0.211324865405187;
00636 
00637         G[m][2][0] = 0.16666666666666666;
00638         G[m][2][1] = 0.16666666666666666;
00639         G[m][2][2] = 0.211324865405187;
00640 
00641         G[m][3][0] = 0.66666666666666666;
00642         G[m][3][1] = 0.16666666666666666;
00643         G[m][3][2] = 0.788675134594813;
00644           
00645         G[m][4][0] = 0.16666666666666666;
00646         G[m][4][1] = 0.66666666666666666;
00647         G[m][4][2] = 0.788675134594813;
00648 
00649         G[m][5][0] = 0.16666666666666666;
00650         G[m][5][1] = 0.16666666666666666;
00651         G[m][5][2] = 0.788675134594813;
00652 
00653         W[m][0] = 0.16666666666666666 / 2.0;
00654         W[m][1] = 0.16666666666666666 / 2.0;
00655         W[m][2] = 0.16666666666666666 / 2.0;
00656         W[m][3] = 0.16666666666666666 / 2.0;
00657         W[m][4] = 0.16666666666666666 / 2.0;
00658         W[m][5] = 0.16666666666666666 / 2.0;
00659           
00660         O[m] = 2;// verify ????????
00661          
00662       }
00663 
00665       FieldVector<double, 3> point(int m, int i)
00666       {
00667         return G[m][i];
00668       }
00669 
00671         double weight (int m, int i)
00672       {
00673         return W[m][i];
00674       }
00675 
00677         int order (int m)
00678       {
00679         return O[m];
00680       }
00681 
00682   private:
00683     FieldVector<double, 3> G[MAXP+1][MAXP];//positions
00684     
00685         double W[MAXP+1][MAXP]; // weights associated with points       
00686         int    O[MAXP+1];       // order of the rule
00687   };
00688 
00689 
00693   template<int dim>  
00694   struct PrismQuadraturePointsSingleton {
00695         static PrismQuadraturePoints<3> prqp;
00696   };
00697 
00701   template<>  
00702   struct PrismQuadraturePointsSingleton<3> {
00703         static PrismQuadraturePoints<3> prqp;
00704   };
00705 
00709   template<typename ct, int dim>
00710   class PrismQuadratureRule;
00711 
00715   template<typename ct>
00716   class PrismQuadratureRule<ct,3> : public QuadratureRule<ct,3>
00717   {
00718   public:
00719 
00721     enum{ d=3 };
00722 
00724     enum{
00725       /* min(Line::order, Triangle::order) */
00726       highest_order =
00727         (int)CubeQuadratureRule<ct,1>::highest_order
00728         < (int)SimplexQuadratureRule<ct,2>::highest_order
00729         ? (int)CubeQuadratureRule<ct,1>::highest_order
00730         : (int)SimplexQuadratureRule<ct,2>::highest_order
00731         };
00732 
00734     typedef ct CoordType;
00735 
00737     typedef PrismQuadratureRule<ct,3> value_type;
00738 
00739     ~PrismQuadratureRule(){}
00740   private:
00741     friend class QuadratureRuleFactory<ct,d>;
00742     PrismQuadratureRule(int p) : QuadratureRule<ct,3>(GeometryType(GeometryType::prism, d))
00743       {
00744         if (p>highest_order)
00745           DUNE_THROW(QuadratureOrderOutOfRange,
00746                      "QuadratureRule for order " << p << " and GeometryType "
00747                      << this->type() << " not available");
00748           
00749         if (p<=2) {
00750           int m=6;
00751           this->delivered_order = PrismQuadraturePointsSingleton<3>::prqp.order(m);
00752           for(int i=0;i<m;++i)
00753           { 
00754             FieldVector<ct,3> local;
00755             for (int k=0; k<d; k++)
00756               local[k] = PrismQuadraturePointsSingleton<3>::prqp.point(m,i)[k];
00757             double weight =
00758               PrismQuadraturePointsSingleton<3>::prqp.weight(m,i);
00759             // put in container
00760             this->push_back(QuadraturePoint<ct,d>(local,weight));
00761           }
00762         }
00763         else {
00764           const QuadratureRule<ct,2> & triangle = QuadratureRules<ct,2>::rule(GeometryType::simplex, p);
00765           const QuadratureRule<ct,1> & line = QuadratureRules<ct,1>::rule(GeometryType::cube, p);
00766           
00767           this->delivered_order = std::min(triangle.order(),line.order());
00768 
00769           for (typename QuadratureRule<ct,1>::const_iterator
00770                  lit = line.begin(); lit != line.end(); ++lit)
00771           {
00772             for (typename QuadratureRule<ct,2>::const_iterator
00773                    tit = triangle.begin(); tit != triangle.end(); ++tit)
00774             {
00775               FieldVector<ct, d> local;
00776               local[0] = tit->position()[0];
00777               local[1] = tit->position()[1];
00778               local[2] = lit->position()[0];
00779               
00780               double weight = tit->weight() * lit->weight();
00781 
00782               // put in container
00783               this->push_back(QuadraturePoint<ct,d>(local,weight));
00784             }
00785           }
00786         }
00787       }
00788   };
00789 
00791   class PyramidQuadraturePoints
00792   {
00793   public:
00794         enum { MAXP=8};
00795         enum { highest_order=2 };
00796 
00798         PyramidQuadraturePoints()
00799       {
00800         int m = 0;
00801         O[m] = 0;
00802          
00803 
00804         // polynom degree 2  ???
00805         m = 8;
00806         G[m][0][0] =0.58541020;
00807         G[m][0][1] =0.72819660;
00808         G[m][0][2] =0.13819660;
00809 
00810         G[m][1][0] =0.13819660;
00811         G[m][1][1] =0.72819660;
00812         G[m][1][2] =0.13819660;
00813 
00814         G[m][2][0] =0.13819660;
00815         G[m][2][1] =0.27630920;
00816         G[m][2][2] =0.58541020;
00817 
00818         G[m][3][0] =0.13819660;
00819         G[m][3][1] =0.27630920;
00820         G[m][3][2] =0.13819660;
00821           
00822         G[m][4][0] =0.72819660;
00823         G[m][4][1] =0.13819660;
00824         G[m][4][2] =0.13819660;
00825 
00826         G[m][5][0] =0.72819660;
00827         G[m][5][1] =0.58541020;
00828         G[m][5][2] =0.13819660;
00829 
00830         G[m][6][0] =0.27630920;
00831         G[m][6][1] =0.13819660;
00832         G[m][6][2] =0.58541020;
00833 
00834         G[m][7][0] =0.27630920;
00835         G[m][7][1] =0.13819660;
00836         G[m][7][2] =0.13819660;
00837 
00838         W[m][0] = 0.125 / 3.0;
00839         W[m][1] = 0.125 / 3.0;
00840         W[m][2] = 0.125 / 3.0;
00841         W[m][3] = 0.125 / 3.0;
00842         W[m][4] = 0.125 / 3.0;
00843         W[m][5] = 0.125 / 3.0;
00844         W[m][6] = 0.125 / 3.0;
00845         W[m][7] = 0.125 / 3.0;
00846           
00847         O[m] = 2;// verify ????????
00848          
00849       }
00850 
00852     FieldVector<double, 3> point(int m, int i)
00853       {
00854         return G[m][i];
00855       }
00856 
00858         double weight (int m, int i)
00859       {
00860         return W[m][i];
00861       }
00862 
00864         int order (int m)
00865       {
00866         return O[m];
00867       }
00868 
00869   private:
00870     FieldVector<double, 3> G[MAXP+1][MAXP];
00871         double W[MAXP+1][MAXP]; // weights associated with points       
00872         int    O[MAXP+1];       // order of the rule
00873   };
00874 
00878   template<int dim>
00879   struct PyramidQuadraturePointsSingleton {};
00880 
00884   template<>
00885   struct PyramidQuadraturePointsSingleton<3> {
00886         static PyramidQuadraturePoints pyqp;
00887   };
00888 
00892   template<typename ct, int dim>
00893   class PyramidQuadratureRule; 
00894 
00898   template<typename ct>
00899   class PyramidQuadratureRule<ct,3> : public QuadratureRule<ct,3>
00900   {
00901   public:
00902 
00904     enum{d=3};
00905 
00907     enum{highest_order=CubeQuadratureRule<ct,1>::highest_order};
00908 
00910     typedef ct CoordType;
00911       
00913     typedef PyramidQuadratureRule<ct,3> value_type;
00914 
00915     ~PyramidQuadratureRule(){}
00916   private:
00917     friend class QuadratureRuleFactory<ct,d>;
00918     PyramidQuadratureRule(int p) : QuadratureRule<ct,3>(GeometryType(GeometryType::pyramid, d))
00919       {
00920         int m;
00921         
00922         if (p>highest_order)
00923           DUNE_THROW(QuadratureOrderOutOfRange,
00924                      "QuadratureRule for order " << p << " and GeometryType "
00925                      << this->type() << " not available");
00926           
00927         if(false) {
00928 //        if(p<=2) {
00929           m=8;
00930           this->delivered_order = PyramidQuadraturePointsSingleton<3>::pyqp.order(m);
00931           FieldVector<ct, d> local;
00932           double weight;
00933           for(int i=0;i<m;++i)
00934           { 
00935             for(int k=0;k<d;++k)
00936               local[k]=PyramidQuadraturePointsSingleton<3>::pyqp.point(m,i)[k];
00937             weight=PyramidQuadraturePointsSingleton<3>::pyqp.weight(m,i);
00938             // put in container
00939             this->push_back(QuadraturePoint<ct,d>(local,weight));       
00940           }
00941         }
00942         else
00943         {
00944           // Define the quadrature rules...
00945           QuadratureRule<ct,3> simplex =
00946             QuadratureRules<ct,3>::rule(GeometryType::simplex,p);
00947 
00948           for (typename QuadratureRule<ct,3>::const_iterator
00949                  it=simplex.begin(); it != simplex.end(); ++it)
00950           {
00951             FieldVector<ct,3> local = it->position();
00952             ct weight = it->weight();
00953             // Simplex 1
00954             // x:=x+y
00955             local[0] = local[0]+local[1];
00956             this->push_back(QuadraturePoint<ct,d>(local,weight));
00957             // Simplex 2
00958             // y:=x+y
00959             local[0] = it->position()[0];
00960             local[1] = local[0]+local[1];
00961             this->push_back(QuadraturePoint<ct,d>(local,weight));
00962           }
00963 
00964           this->delivered_order = simplex.order();
00965         }
00966       }
00967   };
00968 
00975   template<typename ctype, int dim>
00976   class QuadratureRuleFactory {
00977   private:
00978     friend class QuadratureRules<ctype, dim>;
00979     static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
00980       {
00981         if (t.isCube())
00982         {
00983           return CubeQuadratureRule<ctype,dim>(p);
00984         }
00985         if (t.isSimplex())
00986         {
00987           return SimplexQuadratureRule<ctype,dim>(p);
00988         }
00989         DUNE_THROW(Exception, "Unknown GeometryType");
00990       }
00991   };
00992 
00993   template<typename ctype>
00994   class QuadratureRuleFactory<ctype, 0> {
00995   private:
00996     enum { dim = 0 };
00997     friend class QuadratureRules<ctype, dim>;
00998     static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
00999       {
01000         if (t.isVertex())
01001         {
01002           return CubeQuadratureRule<ctype,dim>(p);
01003         }
01004         DUNE_THROW(Exception, "Unknown GeometryType");
01005       }
01006   };
01007 
01008   template<typename ctype>
01009   class QuadratureRuleFactory<ctype, 1> {
01010   private:
01011     enum { dim = 1 };
01012     friend class QuadratureRules<ctype, dim>;
01013     static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
01014       {
01015         if (t.isLine())
01016         {
01017           switch (qt) {
01018           case QuadratureType::Gauss:
01019             return CubeQuadratureRule<ctype,dim>(p);
01020           case QuadratureType::Jacobian_1_0:
01021             return Jacobi1QuadratureRule<ctype,dim>(p);
01022           case QuadratureType::Jacobian_2_0:
01023             return Jacobi2QuadratureRule<ctype,dim>(p);
01024           default:
01025             DUNE_THROW(Exception, "Unknown QuadratureType");
01026           }
01027         }
01028         DUNE_THROW(Exception, "Unknown GeometryType");
01029       }
01030   };
01031 
01032   template<typename ctype>
01033   class QuadratureRuleFactory<ctype, 3> {
01034   private:
01035     enum { dim = 3 };
01036     friend class QuadratureRules<ctype, dim>;
01037     static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
01038       {
01039         if (t.isCube())
01040         {
01041           return CubeQuadratureRule<ctype,dim>(p);
01042         }
01043         if (t.isSimplex())
01044         {
01045           return SimplexQuadratureRule<ctype,dim>(p);
01046         }
01047         if (t.isPrism())
01048             {
01049           return PrismQuadratureRule<ctype,dim>(p);
01050             }
01051         if (t.isPyramid())
01052             {
01053           return PyramidQuadratureRule<ctype,dim>(p);
01054             }
01055         DUNE_THROW(Exception, "Unknown GeometryType");
01056       }
01057   };
01058 
01059 } // end namespace
01060 
01061 #endif