quadraturerules.hh

Go to the documentation of this file.
00001 #ifndef DUNE_QUADRATURERULES_HH
00002 #define DUNE_QUADRATURERULES_HH
00003 
00004 #include<iostream>
00005 #include<vector>
00006 #include<map>
00007 
00008 #include<dune/common/fvector.hh>
00009 #include<dune/common/exceptions.hh>
00010 #include<dune/common/stdstreams.hh>
00011 #include<dune/common/geometrytype.hh>
00012 #include<dune/grid/common/referenceelements.hh>
00013 
00019 namespace Dune {
00020 
00025   class QuadratureOrderOutOfRange : public NotImplemented {};
00026   
00030   template<typename ct, int dim>
00031   class QuadraturePoint {
00032   public:
00033         // compile time parameters
00034         enum { d=dim };
00035         typedef ct CoordType;
00036 
00038         QuadraturePoint (const FieldVector<ct, dim>& x, double w) : local(x)
00039       {
00040         wght = w;
00041       }
00042 
00044         const FieldVector<ct, dim>& position () const
00045       {
00046         return local;
00047       }
00048 
00050         double weight () const
00051       {
00052         return wght;
00053       }
00054     virtual ~QuadraturePoint(){}
00055     
00056   protected:
00057         FieldVector<ct, dim> local;
00058         double wght;
00059   };
00060 
00064   namespace QuadratureType {
00065     enum Enum {
00066       Gauss      = 0,
00067       
00068       Jacobian_1_0 = 1,
00069       Jacobian_2_0 = 2,
00070       
00071       Simpson    = 3,
00072       Trap       = 4,
00073       Grid       = 5,
00074       
00075       Clough     = 21,
00076       
00077       Invalid_Rule = 127
00078     };
00079   };
00080 
00084   template<typename ct, int dim>
00085   class QuadratureRule : public std::vector<QuadraturePoint<ct,dim> >
00086   {
00087   public:
00088 
00090     QuadratureRule() : delivered_order(-1) {}
00091 
00093     QuadratureRule(GeometryType t) : geometry_type(t), delivered_order(-1) {}
00094 
00096     QuadratureRule(GeometryType t, int order) : geometry_type(t), delivered_order(order) {}
00097 
00099         enum { d=dim };
00100 
00102         typedef ct CoordType;
00103 
00105         virtual int order () const { return delivered_order; }
00106 
00108         virtual GeometryType type () const { return geometry_type; }
00109     virtual ~QuadratureRule(){}
00110 
00113     typedef typename std::vector<QuadraturePoint<ct,dim> >::const_iterator iterator;
00114 
00115   protected:
00116     GeometryType geometry_type;
00117         int delivered_order;
00118     
00119     void tensor_product (const QuadratureRule<ct,1> & q1d)
00120       {
00121         // fill in all the gauss points
00122         int m = q1d.size();
00123         int n = power(m,dim);
00124         for (int i=0; i<n; i++)
00125                 {
00126                   // compute multi index for Gauss point
00127                   int x[dim];
00128                   int z = i;
00129                   for (int k=0; k<dim; ++k)
00130           {
00131             x[k] = z%m;
00132             z = z/m;
00133           }
00134 
00135                   // compute coordinates and weight
00136                   double weight = 1.0;
00137                   FieldVector<ct, dim> local;
00138                   for (int k=0; k<dim; k++) 
00139           {
00140             local[k] = q1d[x[k]].position()[0];
00141             weight *= q1d[x[k]].weight();
00142           }
00143 
00144                   // put in container
00145                   push_back(QuadraturePoint<ct,dim>(local,weight));
00146                 }
00147       }
00148 
00149         int power (int y, int d)
00150       {
00151         int m=1;
00152         for (int i=0; i<d; i++) m *= y;
00153         return m;
00154       }
00155   };
00156 
00157   // Forward declaration of the factory class,
00158   // needed internally by the QuadratureRules container class.
00159   template<typename ctype, int dim> struct QuadratureRuleFactory;
00160 
00164   template<typename ctype, int dim>
00165   class QuadratureRules {
00166 
00168     typedef std::pair<GeometryType,int> QuadratureRuleKey;
00169 
00171     typedef QuadratureRule<ctype, dim> QuadratureRule;
00172 
00174     const QuadratureRule& _rule(const GeometryType& t, int p, QuadratureType::Enum qt=QuadratureType::Gauss)
00175       {
00176         static std::map<QuadratureRuleKey, QuadratureRule> _quadratureMap;
00177         QuadratureRuleKey key(t,p);
00178         if (_quadratureMap.find(key) == _quadratureMap.end()) {
00179           /*
00180             The rule must be acquired before we can store it.
00181             If we write this in one command, an invalid rule
00182             would get stored in case of an exception.
00183           */
00184           QuadratureRule rule =
00185             QuadratureRuleFactory<ctype,dim>::rule(t,p,qt);
00186           _quadratureMap[key] = rule;
00187         }
00188         return _quadratureMap[key];
00189       }
00191     static QuadratureRules& instance()
00192       {
00193         static QuadratureRules instance;
00194         return instance;
00195       }
00197     QuadratureRules () {};
00198   public:
00200     static const QuadratureRule& rule(const GeometryType& t, int p, QuadratureType::Enum qt=QuadratureType::Gauss)
00201       {
00202         return instance()._rule(t,p,qt);
00203       }
00205     static const QuadratureRule& rule(const GeometryType::BasicType t, int p, QuadratureType::Enum qt=QuadratureType::Gauss)
00206       {
00207         GeometryType gt(t,dim);
00208         return instance()._rule(gt,p,qt);
00209       }
00210   };
00211 
00212   /***********************************************************/
00213 
00217   template<typename ct, int dim>
00218   class CubeQuadratureRule :
00219     public QuadratureRule<ct,dim>
00220   {
00221   public:
00223         enum { d=dim };
00224 
00226         enum { highest_order=CubeQuadratureRule<ct,1>::highest_order };
00227 
00229         typedef ct CoordType;
00230 
00232         typedef CubeQuadratureRule value_type;
00233 
00234     ~CubeQuadratureRule(){}
00235   private:
00236     friend class QuadratureRuleFactory<ct,dim>;
00238         CubeQuadratureRule (int p) : QuadratureRule<ct,dim>(GeometryType(GeometryType::cube, d))
00239       {
00240         QuadratureRule<ct,1> q1D = QuadratureRules<ct,1>::rule(GeometryType::cube, p);
00241         tensor_product( q1D );
00242         this->delivered_order = q1D.order();
00243       }
00244 
00245   };
00246 
00249   template<typename ct>
00250   class CubeQuadratureRule<ct,0> :
00251     public QuadratureRule<ct,0>
00252   {
00253   public:
00254         // compile time parameters
00255         enum { d=0 };
00256         enum { dim=0 };
00257         enum { highest_order=1000 };
00258         typedef ct CoordType;
00259         typedef CubeQuadratureRule value_type;
00260 
00261     ~CubeQuadratureRule(){}
00262   private:
00263     friend class QuadratureRuleFactory<ct,dim>;
00264     CubeQuadratureRule (int p):
00265       QuadratureRule<ct,0>(GeometryType(GeometryType::cube, 0))
00266     {
00267       FieldVector<ct, dim> point(0.0);
00268 
00269       if (p > highest_order)
00270         DUNE_THROW(QuadratureOrderOutOfRange, "Quadrature rule " << p << " not supported!");
00271 
00272       this->delivered_order = highest_order;
00273       this->push_back(QuadraturePoint<ct,dim>(point, 1.0));
00274     }
00275   };
00276  
00279   template<typename ct>
00280   class CubeQuadratureRule<ct,1> :
00281     public QuadratureRule<ct,1>
00282   {
00283   public:
00284         // compile time parameters
00285         enum { d=1 };
00286         enum { dim=1 };
00287         enum { highest_order=44 };
00288         typedef ct CoordType;
00289         typedef CubeQuadratureRule value_type;
00290 
00291     ~CubeQuadratureRule(){}
00292   private:
00293     void init(int p,
00294               std::vector< FieldVector<ct, dim> > & _points,
00295               std::vector< double > & _weight,
00296               int & delivered_order);
00297     friend class QuadratureRuleFactory<ct,dim>;
00298     CubeQuadratureRule (int p)
00299       : QuadratureRule<ct,1>(GeometryType(GeometryType::cube, 1))
00300     {
00302       std::vector< FieldVector<ct, dim> > _points;
00303       std::vector< double > _weight;
00304       
00305       int delivered_order;
00306       
00307       init(p, _points, _weight, delivered_order);
00308       
00309       this->delivered_order = delivered_order;
00310       assert(_points.size() == _weight.size());
00311       for (size_t i = 0; i < _points.size(); i++)
00312         this->push_back(QuadraturePoint<ct,dim>(_points[i], _weight[i]));
00313     }
00314   };
00315  
00319   template<typename ct, int dim>
00320   class Jacobi1QuadratureRule;
00321   
00325   template<typename ct>
00326   class Jacobi1QuadratureRule<ct,1> :
00327     public QuadratureRule<ct,1>
00328   {
00329   public:
00331         enum { d=1 };
00334         enum { dim=1 };
00335 
00337         enum { highest_order=44 };
00338 
00340         typedef ct CoordType;
00341 
00343         typedef Jacobi1QuadratureRule value_type;
00344 
00345     ~Jacobi1QuadratureRule(){}
00346   private:
00347     void init(int p,
00348               std::vector< FieldVector<ct, dim> > & _points,
00349               std::vector< double > & _weight,
00350               int & delivered_order);
00351     friend class QuadratureRuleFactory<ct,dim>;
00352     Jacobi1QuadratureRule (int p)
00353       : QuadratureRule<ct,1>(GeometryType(GeometryType::cube, 1))
00354     {
00356       std::vector< FieldVector<ct, dim> > _points;
00357       std::vector< double > _weight;
00358       
00359       int delivered_order;
00360       
00361       init(p, _points, _weight, delivered_order);
00362       this->delivered_order = delivered_order;
00363       assert(_points.size() == _weight.size());
00364       for (size_t i = 0; i < _points.size(); i++)
00365         this->push_back(QuadraturePoint<ct,dim>(_points[i], _weight[i]));
00366     }
00367   };
00368  
00372   template<typename ct, int dim>
00373   class Jacobi2QuadratureRule;
00374   
00378   template<typename ct>
00379   class Jacobi2QuadratureRule<ct,1> :
00380     public QuadratureRule<ct,1>
00381   {
00382   public:
00384         enum { d=1 };
00385 
00388         enum { dim=1 };
00389 
00391         enum { highest_order=44 };
00392 
00394         typedef ct CoordType;
00395 
00397         typedef Jacobi2QuadratureRule value_type;
00398 
00399     ~Jacobi2QuadratureRule(){}
00400   private:
00401     void init(int p,
00402               std::vector< FieldVector<ct, dim> > & _points,
00403               std::vector< double > & _weight,
00404               int & delivered_order);
00405     friend class QuadratureRuleFactory<ct,dim>;
00406     Jacobi2QuadratureRule (int p)
00407       : QuadratureRule<ct,1>(GeometryType(GeometryType::cube, 1))
00408     {
00410       std::vector< FieldVector<ct, dim> > _points;
00411       std::vector< double > _weight;
00412       
00413       int delivered_order;
00414       
00415       init(p, _points, _weight, delivered_order);
00416       
00417       this->delivered_order = delivered_order;
00418       assert(_points.size() == _weight.size());
00419       for (size_t i = 0; i < _points.size(); i++)
00420         this->push_back(QuadraturePoint<ct,dim>(_points[i], _weight[i]));
00421     }
00422   };
00423  
00424   /************************************************
00425    * Quadraturerule for Simplices/Triangle
00426    *************************************************/
00427 
00431   template<typename ct, int dim>
00432   class SimplexQuadratureRule;
00433 
00437   template<typename ct>
00438   class SimplexQuadratureRule<ct,2> : public QuadratureRule<ct,2>
00439   {
00440   public:
00441 
00443     enum{d=2};
00444 
00446     enum { highest_order=44 };
00447 
00449     typedef ct CoordType;
00450 
00452     typedef SimplexQuadratureRule value_type;
00453     ~SimplexQuadratureRule(){}
00454   private:
00455     friend class QuadratureRuleFactory<ct,d>;
00456     SimplexQuadratureRule (int p);
00457   };
00458 
00462   template<typename ct>
00463   class SimplexQuadratureRule<ct,3> : public QuadratureRule<ct,3>
00464   {
00465   public:
00466 
00468     enum{d=3};
00469 
00471     enum { highest_order=44 };
00472 
00474     typedef ct CoordType;
00475 
00477     typedef SimplexQuadratureRule<ct,3> value_type;
00478     ~SimplexQuadratureRule(){}
00479   private:
00480     friend class QuadratureRuleFactory<ct,d>;
00481     SimplexQuadratureRule (int p);
00482   };
00483 
00484 /***********************************
00485  * quadrature for Prism
00486  **********************************/
00487 
00489   template<int dim>
00490   class PrismQuadraturePoints;
00491 
00493   template<>
00494   class PrismQuadraturePoints<3>
00495   {
00496   public:
00497         enum { MAXP=6};
00498         enum { highest_order=2 };
00499 
00501         PrismQuadraturePoints ()
00502       {
00503         int m = 0;
00504         O[m] = 0;
00505           
00506         // polynom degree 0  ???
00507         m = 6;
00508         G[m][0][0] = 0.0;
00509         G[m][0][1] = 0.0;
00510         G[m][0][2] = 0.0;
00511 
00512         G[m][1][0] = 1.0;
00513         G[m][1][1] = 0.0;
00514         G[m][1][2] = 0.0;
00515 
00516         G[m][2][0] = 0.0;
00517         G[m][2][1] = 1.0;
00518         G[m][2][2] = 0.0;
00519 
00520         G[m][3][0] = 0.0;
00521         G[m][3][1] = 0.0;
00522         G[m][3][2] = 1.0;
00523           
00524         G[m][4][0] = 1.0;
00525         G[m][4][1] = 0.0;
00526         G[m][4][2] = 1.0;
00527 
00528         G[m][5][0] = 0.0;
00529         G[m][5][1] = 0.1;
00530         G[m][5][2] = 1.0;
00531 
00532         W[m][0] = 0.16666666666666666 / 2.0;
00533         W[m][1] = 0.16666666666666666 / 2.0;
00534         W[m][2] = 0.16666666666666666 / 2.0;
00535         W[m][3] = 0.16666666666666666 / 2.0;
00536         W[m][4] = 0.16666666666666666 / 2.0;
00537         W[m][5] = 0.16666666666666666 / 2.0;
00538           
00539         O[m] = 0;// verify ????????
00540           
00541 
00542         // polynom degree 2  ???
00543         m = 6;
00544         G[m][0][0] =0.66666666666666666 ;
00545         G[m][0][1] =0.16666666666666666 ;
00546         G[m][0][2] =0.211324865405187 ;
00547 
00548         G[m][1][0] = 0.16666666666666666;
00549         G[m][1][1] =0.66666666666666666 ;
00550         G[m][1][2] = 0.211324865405187;
00551 
00552         G[m][2][0] = 0.16666666666666666;
00553         G[m][2][1] = 0.16666666666666666;
00554         G[m][2][2] = 0.211324865405187;
00555 
00556         G[m][3][0] = 0.66666666666666666;
00557         G[m][3][1] = 0.16666666666666666;
00558         G[m][3][2] = 0.788675134594813;
00559           
00560         G[m][4][0] = 0.16666666666666666;
00561         G[m][4][1] = 0.66666666666666666;
00562         G[m][4][2] = 0.788675134594813;
00563 
00564         G[m][5][0] = 0.16666666666666666;
00565         G[m][5][1] = 0.16666666666666666;
00566         G[m][5][2] = 0.788675134594813;
00567 
00568         W[m][0] = 0.16666666666666666 / 2.0;
00569         W[m][1] = 0.16666666666666666 / 2.0;
00570         W[m][2] = 0.16666666666666666 / 2.0;
00571         W[m][3] = 0.16666666666666666 / 2.0;
00572         W[m][4] = 0.16666666666666666 / 2.0;
00573         W[m][5] = 0.16666666666666666 / 2.0;
00574           
00575         O[m] = 2;// verify ????????
00576          
00577       }
00578 
00580       FieldVector<double, 3> point(int m, int i)
00581       {
00582         return G[m][i];
00583       }
00584 
00586         double weight (int m, int i)
00587       {
00588         return W[m][i];
00589       }
00590 
00592         int order (int m)
00593       {
00594         return O[m];
00595       }
00596 
00597   private:
00598     FieldVector<double, 3> G[MAXP+1][MAXP];//positions
00599     
00600         double W[MAXP+1][MAXP]; // weights associated with points       
00601         int    O[MAXP+1];       // order of the rule
00602   };
00603 
00604 
00608   template<int dim>  
00609   struct PrismQuadraturePointsSingleton {
00610         static PrismQuadraturePoints<3> prqp;
00611   };
00612 
00616   template<>  
00617   struct PrismQuadraturePointsSingleton<3> {
00618         static PrismQuadraturePoints<3> prqp;
00619   };
00620 
00624   template<typename ct, int dim>
00625   class PrismQuadratureRule;
00626 
00630   template<typename ct>
00631   class PrismQuadratureRule<ct,3> : public QuadratureRule<ct,3>
00632   {
00633   public:
00634 
00636     enum{ d=3 };
00637 
00639     enum{
00640       /* min(Line::order, Triangle::order) */
00641       highest_order =
00642         (int)CubeQuadratureRule<ct,1>::highest_order
00643         < (int)SimplexQuadratureRule<ct,2>::highest_order
00644         ? (int)CubeQuadratureRule<ct,1>::highest_order
00645         : (int)SimplexQuadratureRule<ct,2>::highest_order
00646         };
00647 
00649     typedef ct CoordType;
00650 
00652     typedef PrismQuadratureRule<ct,3> value_type;
00653 
00654     ~PrismQuadratureRule(){}
00655   private:
00656     friend class QuadratureRuleFactory<ct,d>;
00657     PrismQuadratureRule(int p) : QuadratureRule<ct,3>(GeometryType(GeometryType::prism, d))
00658       {
00659         if (p>highest_order)
00660           DUNE_THROW(QuadratureOrderOutOfRange,
00661                      "QuadratureRule for order " << p << " and GeometryType "
00662                      << this->type() << " not available");
00663           
00664         if (p<=2) {
00665           int m=6;
00666           this->delivered_order = PrismQuadraturePointsSingleton<3>::prqp.order(m);
00667           for(int i=0;i<m;++i)
00668           { 
00669             FieldVector<ct,3> local;
00670             for (int k=0; k<d; k++)
00671               local[k] = PrismQuadraturePointsSingleton<3>::prqp.point(m,i)[k];
00672             double weight =
00673               PrismQuadraturePointsSingleton<3>::prqp.weight(m,i);
00674             // put in container
00675             push_back(QuadraturePoint<ct,d>(local,weight));
00676           }
00677         }
00678         else {
00679           const QuadratureRule<ct,2> & triangle = QuadratureRules<ct,2>::rule(GeometryType::simplex, p);
00680           const QuadratureRule<ct,1> & line = QuadratureRules<ct,1>::rule(GeometryType::cube, p);
00681           
00682           this->delivered_order = std::min(triangle.order(),line.order());
00683 
00684           for (typename QuadratureRule<ct,1>::const_iterator
00685                  lit = line.begin(); lit != line.end(); ++lit)
00686           {
00687             for (typename QuadratureRule<ct,2>::const_iterator
00688                    tit = triangle.begin(); tit != triangle.end(); ++tit)
00689             {
00690               FieldVector<ct, d> local;
00691               local[0] = tit->position()[0];
00692               local[1] = tit->position()[1];
00693               local[2] = lit->position()[0];
00694               
00695               double weight = tit->weight() * lit->weight();
00696 
00697               // put in container
00698               push_back(QuadraturePoint<ct,d>(local,weight));
00699             }
00700           }
00701         }
00702       }
00703   };
00704 
00706   class PyramidQuadraturePoints
00707   {
00708   public:
00709         enum { MAXP=8};
00710         enum { highest_order=2 };
00711 
00713         PyramidQuadraturePoints()
00714       {
00715         int m = 0;
00716         O[m] = 0;
00717          
00718 
00719         // polynom degree 2  ???
00720         m = 8;
00721         G[m][0][0] =0.58541020;
00722         G[m][0][1] =0.72819660;
00723         G[m][0][2] =0.13819660;
00724 
00725         G[m][1][0] =0.13819660;
00726         G[m][1][1] =0.72819660;
00727         G[m][1][2] =0.13819660;
00728 
00729         G[m][2][0] =0.13819660;
00730         G[m][2][1] =0.27630920;
00731         G[m][2][2] =0.58541020;
00732 
00733         G[m][3][0] =0.13819660;
00734         G[m][3][1] =0.27630920;
00735         G[m][3][2] =0.13819660;
00736           
00737         G[m][4][0] =0.72819660;
00738         G[m][4][1] =0.13819660;
00739         G[m][4][2] =0.13819660;
00740 
00741         G[m][5][0] =0.72819660;
00742         G[m][5][1] =0.58541020;
00743         G[m][5][2] =0.13819660;
00744 
00745         G[m][6][0] =0.27630920;
00746         G[m][6][1] =0.13819660;
00747         G[m][6][2] =0.58541020;
00748 
00749         G[m][7][0] =0.27630920;
00750         G[m][7][1] =0.13819660;
00751         G[m][7][2] =0.13819660;
00752 
00753         W[m][0] = 0.125 / 3.0;
00754         W[m][1] = 0.125 / 3.0;
00755         W[m][2] = 0.125 / 3.0;
00756         W[m][3] = 0.125 / 3.0;
00757         W[m][4] = 0.125 / 3.0;
00758         W[m][5] = 0.125 / 3.0;
00759         W[m][6] = 0.125 / 3.0;
00760         W[m][7] = 0.125 / 3.0;
00761           
00762         O[m] = 2;// verify ????????
00763          
00764       }
00765 
00767     FieldVector<double, 3> point(int m, int i)
00768       {
00769         return G[m][i];
00770       }
00771 
00773         double weight (int m, int i)
00774       {
00775         return W[m][i];
00776       }
00777 
00779         int order (int m)
00780       {
00781         return O[m];
00782       }
00783 
00784   private:
00785     FieldVector<double, 3> G[MAXP+1][MAXP];
00786         double W[MAXP+1][MAXP]; // weights associated with points       
00787         int    O[MAXP+1];       // order of the rule
00788   };
00789 
00793   template<int dim>
00794   struct PyramidQuadraturePointsSingleton {};
00795 
00799   template<>
00800   struct PyramidQuadraturePointsSingleton<3> {
00801         static PyramidQuadraturePoints pyqp;
00802   };
00803 
00807   template<typename ct, int dim>
00808   class PyramidQuadratureRule; 
00809 
00813   template<typename ct>
00814   class PyramidQuadratureRule<ct,3> : public QuadratureRule<ct,3>
00815   {
00816   public:
00817 
00819     enum{d=3};
00820 
00822     enum{highest_order=CubeQuadratureRule<ct,1>::highest_order};
00823 
00825     typedef ct CoordType;
00826       
00828     typedef PyramidQuadratureRule<ct,3> value_type;
00829 
00830     ~PyramidQuadratureRule(){}
00831   private:
00832     friend class QuadratureRuleFactory<ct,d>;
00833     PyramidQuadratureRule(int p) : QuadratureRule<ct,3>(GeometryType(GeometryType::pyramid, d))
00834       {
00835         int m;
00836         
00837         if (p>highest_order)
00838           DUNE_THROW(QuadratureOrderOutOfRange,
00839                      "QuadratureRule for order " << p << " and GeometryType "
00840                      << this->type() << " not available");
00841           
00842         if(false) {
00843 //        if(p<=2) {
00844           m=8;
00845           this->delivered_order = PyramidQuadraturePointsSingleton<3>::pyqp.order(m);
00846           FieldVector<ct, d> local;
00847           double weight;
00848           for(int i=0;i<m;++i)
00849           { 
00850             for(int k=0;k<d;++k)
00851               local[k]=PyramidQuadraturePointsSingleton<3>::pyqp.point(m,i)[k];
00852             weight=PyramidQuadraturePointsSingleton<3>::pyqp.weight(m,i);
00853             // put in container
00854             push_back(QuadraturePoint<ct,d>(local,weight));     
00855           }
00856         }
00857         else
00858         {
00859           // Define the quadrature rules...
00860           QuadratureRule<ct,3> simplex =
00861             QuadratureRules<ct,3>::rule(GeometryType::simplex,p);
00862 
00863           for (typename QuadratureRule<ct,3>::const_iterator
00864                  it=simplex.begin(); it != simplex.end(); ++it)
00865           {
00866             FieldVector<ct,3> local = it->position();
00867             ct weight = it->weight();
00868             // Simplex 1
00869             // x:=x+y
00870             local[0] = local[0]+local[1];
00871             push_back(QuadraturePoint<ct,d>(local,weight));
00872             // Simplex 2
00873             // y:=x+y
00874             local[0] = it->position()[0];
00875             local[1] = local[0]+local[1];
00876             push_back(QuadraturePoint<ct,d>(local,weight));
00877           }
00878 
00879           this->delivered_order = simplex.order();
00880         }
00881       }
00882   };
00883 
00890   template<typename ctype, int dim>
00891   class QuadratureRuleFactory {
00892   private:
00893     friend class QuadratureRules<ctype, dim>;
00894     static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
00895       {
00896         if (t.isCube())
00897         {
00898           return CubeQuadratureRule<ctype,dim>(p);
00899         }
00900         if (t.isSimplex())
00901         {
00902           return SimplexQuadratureRule<ctype,dim>(p);
00903         }
00904         DUNE_THROW(Exception, "Unknown GeometryType");
00905       }
00906   };
00907 
00908   template<typename ctype>
00909   class QuadratureRuleFactory<ctype, 0> {
00910   private:
00911     enum { dim = 0 };
00912     friend class QuadratureRules<ctype, dim>;
00913     static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
00914       {
00915         if (t.isVertex())
00916         {
00917           return CubeQuadratureRule<ctype,dim>(p);
00918         }
00919         DUNE_THROW(Exception, "Unknown GeometryType");
00920       }
00921   };
00922 
00923   template<typename ctype>
00924   class QuadratureRuleFactory<ctype, 1> {
00925   private:
00926     enum { dim = 1 };
00927     friend class QuadratureRules<ctype, dim>;
00928     static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
00929       {
00930         if (t.isLine())
00931         {
00932           switch (qt) {
00933           case QuadratureType::Gauss:
00934             return CubeQuadratureRule<ctype,dim>(p);
00935           case QuadratureType::Jacobian_1_0:
00936             return Jacobi1QuadratureRule<ctype,dim>(p);
00937           case QuadratureType::Jacobian_2_0:
00938             return Jacobi2QuadratureRule<ctype,dim>(p);
00939           default:
00940             DUNE_THROW(Exception, "Unknown QuadratureType");
00941           }
00942         }
00943         DUNE_THROW(Exception, "Unknown GeometryType");
00944       }
00945   };
00946 
00947   template<typename ctype>
00948   class QuadratureRuleFactory<ctype, 3> {
00949   private:
00950     enum { dim = 3 };
00951     friend class QuadratureRules<ctype, dim>;
00952     static QuadratureRule<ctype, dim> rule(const GeometryType& t, int p, QuadratureType::Enum qt)
00953       {
00954         if (t.isCube())
00955         {
00956           return CubeQuadratureRule<ctype,dim>(p);
00957         }
00958         if (t.isSimplex())
00959         {
00960           return SimplexQuadratureRule<ctype,dim>(p);
00961         }
00962         if (t.isPrism())
00963             {
00964           return PrismQuadratureRule<ctype,dim>(p);
00965             }
00966         if (t.isPyramid())
00967             {
00968           return PyramidQuadratureRule<ctype,dim>(p);
00969             }
00970         DUNE_THROW(Exception, "Unknown GeometryType");
00971       }
00972   };
00973 
00974 } // end namespace
00975 
00976 #endif

Generated on 12 Dec 2007 with Doxygen (ver 1.5.1)