dune-localfunctions  2.1.1
monomialbasis.hh
Go to the documentation of this file.
00001 #ifndef DUNE_MONOMIALBASIS_HH
00002 #define DUNE_MONOMIALBASIS_HH
00003 
00004 #include <vector>
00005 
00006 #include <dune/common/fvector.hh>
00007 #include <dune/common/fmatrix.hh>
00008 
00009 #include <dune/grid/genericgeometry/topologytypes.hh>
00010 
00011 #include <dune/localfunctions/utility/field.hh>
00012 
00013 #include <dune/grid/common/topologyfactory.hh>
00014 #include <dune/localfunctions/utility/multiindex.hh>
00015 #include <dune/localfunctions/utility/tensor.hh>
00016 
00017 namespace Dune
00018 {
00019 /************************************************
00020  * Classes for evaluating ''Monomials'' on any order
00021  * for all reference element type. 
00022  * For a simplex topology these are the nomral 
00023  * monomials for cube topologies the bimonomials.
00024  * The construction follows the construction of the
00025  * generic geometries using tensor products for
00026  * prism generation and duffy transform for pyramid
00027  * construction.
00028  * A derivative argument can be applied, in which case
00029  * all derivatives up to the desired order are 
00030  * evaluated. Note that in for higher order derivatives
00031  * only the ''lower'' part of the symmetric tensor
00032  * is evaluated, e.g., passing derivative equal to 2
00033  * to the class will provide the vector
00034  *    (d/dxdx p, d/dxydx p, d/dydy p, 
00035  *     d/dx p, d/dy p, p)
00036  * Important:
00037  * So far the computation of the derivatives has not
00038  * been fully implemented for general pyramid
00039  * construction, i.e., in the case where a pyramid is
00040  * build over a non simplex base geometry.
00041  *
00042  * Central classes:
00043  * 1) template< class Topology, class F >
00044  *    class MonomialBasisImpl;
00045  *    Implementation of the monomial evaluation for
00046  *    a given topology and field type.
00047  *    The method evaluate fills a F* vector
00048  * 2) template< class Topology, class F >
00049  *    class MonomialBasis
00050  *    The base class for the static monomial evaluation
00051  *    providing addiional evaluate methods including
00052  *    one taking std::vector<F>.
00053  * 3) template< int dim, class F >
00054  *    class VirtualMonomialBasis
00055  *    Virtualization of the MonomialBasis.
00056  * 4) template< int dim, class F >
00057  *    struct MonomialBasisFactory;
00058  *    A factory class for the VirtualMonomialBasis
00059  * 5) template< int dim, class F >
00060  *    struct MonomialBasisProvider
00061  *    A singleton container for the virtual monomial
00062  *    basis
00063  ************************************************/
00064 
00065   // Internal Forward Declarations
00066   // -----------------------------
00067   
00068   template< class Topology >
00069   class MonomialBasisSize;
00070 
00071   template< class Topology, class F >
00072   class MonomialBasis;
00073 
00074 
00075 
00076   // MonomialBasisSize
00077   // -----------------
00078 
00079   template<>
00080   class MonomialBasisSize< GenericGeometry::Point >
00081   {
00082     typedef MonomialBasisSize< GenericGeometry::Point > This;
00083 
00084   public:
00085     static This &instance ()
00086     {
00087       static This _instance;
00088       return _instance;
00089     }
00090 
00091     typedef GenericGeometry::Point Topology;
00092 
00093     friend class MonomialBasisSize< GenericGeometry::Prism< Topology > >;
00094     friend class MonomialBasisSize< GenericGeometry::Pyramid< Topology > >;
00095 
00096     mutable unsigned int maxOrder_;
00097     // sizes_[ k ]: number of basis functions of exactly order k
00098     mutable unsigned int *sizes_;
00099     // numBaseFunctions_[ k ] = sizes_[ 0 ] + ... + sizes_[ k ]
00100     mutable unsigned int *numBaseFunctions_;
00101 
00102     MonomialBasisSize ()
00103     : maxOrder_( 0 ),
00104       sizes_( 0 ),
00105       numBaseFunctions_( 0 )
00106     {
00107       computeSizes( 2 );
00108     }
00109 
00110     ~MonomialBasisSize ()
00111     {
00112       delete[] sizes_;
00113       delete[] numBaseFunctions_;
00114     }
00115 
00116     unsigned int operator() ( const unsigned int order ) const
00117     {
00118       return numBaseFunctions_[ order ];
00119     }
00120 
00121     unsigned int maxOrder () const
00122     {
00123       return maxOrder_;
00124     }
00125 
00126     void computeSizes ( unsigned int order ) const
00127     {
00128       if (order <= maxOrder_)
00129         return;
00130 
00131       maxOrder_ = order;
00132 
00133       delete [] sizes_;
00134       delete [] numBaseFunctions_;
00135       sizes_            = new unsigned int [ order+1 ];
00136       numBaseFunctions_ = new unsigned int [ order+1 ];
00137 
00138       sizes_[ 0 ] = 1;
00139       numBaseFunctions_[ 0 ] = 1;
00140       for( unsigned int k = 1; k <= order; ++k )
00141       {
00142         sizes_[ k ]            = 0;
00143         numBaseFunctions_[ k ] = 1;
00144       }
00145     }
00146   };
00147 
00148   template< class BaseTopology >
00149   class MonomialBasisSize< GenericGeometry::Prism< BaseTopology > >
00150   {
00151     typedef MonomialBasisSize< GenericGeometry::Prism< BaseTopology > > This;
00152 
00153   public:
00154     static This &instance ()
00155     {
00156       static This _instance;
00157       return _instance;
00158     }
00159 
00160     typedef GenericGeometry::Prism< BaseTopology > Topology;
00161 
00162     friend class MonomialBasisSize< GenericGeometry::Prism< Topology > >;
00163     friend class MonomialBasisSize< GenericGeometry::Pyramid< Topology > >;
00164 
00165     mutable unsigned int maxOrder_;
00166     // sizes_[ k ]: number of basis functions of exactly order k
00167     mutable unsigned int *sizes_;
00168     // numBaseFunctions_[ k ] = sizes_[ 0 ] + ... + sizes_[ k ]
00169     mutable unsigned int *numBaseFunctions_;
00170 
00171     MonomialBasisSize ()
00172     : maxOrder_( 0 ),
00173       sizes_( 0 ),
00174       numBaseFunctions_( 0 )
00175     {
00176       computeSizes( 2 );
00177     }
00178 
00179     ~MonomialBasisSize ()
00180     {
00181       delete[] sizes_;
00182       delete[] numBaseFunctions_;
00183     }
00184 
00185     unsigned int operator() ( const unsigned int order ) const
00186     {
00187       return numBaseFunctions_[ order ];
00188     }
00189 
00190     unsigned int maxOrder() const
00191     {
00192       return maxOrder_;
00193     }
00194 
00195     void computeSizes ( unsigned int order ) const
00196     {
00197       if (order <= maxOrder_)
00198         return;
00199 
00200       maxOrder_ = order;
00201 
00202       delete[] sizes_;
00203       delete[] numBaseFunctions_;
00204       sizes_            = new unsigned int[ order+1 ];
00205       numBaseFunctions_ = new unsigned int[ order+1 ];
00206 
00207       MonomialBasisSize<BaseTopology> &baseBasis =
00208             MonomialBasisSize<BaseTopology>::instance();
00209       baseBasis.computeSizes( order );
00210       const unsigned int *const baseSizes = baseBasis.sizes_;
00211       const unsigned int *const baseNBF   = baseBasis.numBaseFunctions_;
00212 
00213       sizes_[ 0 ] = 1;
00214       numBaseFunctions_[ 0 ] = 1;
00215       for( unsigned int k = 1; k <= order; ++k )
00216       {
00217         sizes_[ k ]            = baseNBF[ k ] + k*baseSizes[ k ];
00218         numBaseFunctions_[ k ] = numBaseFunctions_[ k-1 ] + sizes_[ k ];
00219       }
00220     }
00221   };
00222 
00223   template< class BaseTopology >
00224   class MonomialBasisSize< GenericGeometry::Pyramid< BaseTopology > >
00225   {
00226     typedef MonomialBasisSize< GenericGeometry::Pyramid< BaseTopology > > This;
00227 
00228   public:
00229     static This &instance ()
00230     {
00231       static This _instance;
00232       return _instance;
00233     }
00234 
00235     typedef GenericGeometry::Pyramid< BaseTopology > Topology;
00236 
00237     friend class MonomialBasisSize< GenericGeometry::Prism< Topology > >;
00238     friend class MonomialBasisSize< GenericGeometry::Pyramid< Topology > >;
00239 
00240     mutable unsigned int maxOrder_;
00241     // sizes_[ k ]: number of basis functions of exactly order k
00242     mutable unsigned int *sizes_;
00243     // numBaseFunctions_[ k ] = sizes_[ 0 ] + ... + sizes_[ k ]
00244     mutable unsigned int *numBaseFunctions_;
00245 
00246     MonomialBasisSize ()
00247     : maxOrder_( 0 ),
00248       sizes_( 0 ),
00249       numBaseFunctions_( 0 )
00250     {
00251       computeSizes( 2 );
00252     }
00253 
00254     ~MonomialBasisSize ()
00255     {
00256       delete[] sizes_;
00257       delete[] numBaseFunctions_;
00258     }
00259 
00260     unsigned int operator() ( const unsigned int order ) const
00261     {
00262       return numBaseFunctions_[ order ];
00263     }
00264 
00265     unsigned int maxOrder() const
00266     {
00267       return maxOrder_;
00268     }
00269 
00270     void computeSizes ( unsigned int order ) const
00271     {
00272       if (order <= maxOrder_)
00273         return;
00274 
00275       maxOrder_ = order;
00276 
00277       delete[] sizes_;
00278       delete[] numBaseFunctions_;
00279       sizes_            = new unsigned int[ order+1 ];
00280       numBaseFunctions_ = new unsigned int[ order+1 ];
00281 
00282       MonomialBasisSize<BaseTopology> &baseBasis =
00283             MonomialBasisSize<BaseTopology>::instance();
00284 
00285       baseBasis.computeSizes( order );
00286 
00287       const unsigned int *const baseNBF = baseBasis.numBaseFunctions_;
00288       sizes_[ 0 ] = 1;
00289       numBaseFunctions_[ 0 ] = 1;
00290       for( unsigned int k = 1; k <= order; ++k )
00291       {
00292         sizes_[ k ]            = baseNBF[ k ];
00293         numBaseFunctions_[ k ] = numBaseFunctions_[ k-1 ] + sizes_[ k ];
00294       }
00295     }
00296   };
00297 
00298 
00299 
00300   // MonomialBasisHelper
00301   // -------------------
00302 
00303 
00304   template< int mydim, int dim, class F >
00305   struct MonomialBasisHelper
00306   {
00307     typedef MonomialBasisSize< typename GenericGeometry::SimplexTopology< mydim >::type > MySize;
00308     typedef MonomialBasisSize< typename GenericGeometry::SimplexTopology< dim >::type > Size;
00309 
00310     static void copy ( const unsigned int deriv, F *&wit, F *&rit,
00311                        const unsigned int numBaseFunctions, const F &z )
00312     {
00313       // n(d,k) = size<k>[d];
00314       MySize &mySize = MySize::instance();
00315       Size &size = Size::instance();
00316 
00317       const F *const rend = rit + size( deriv )*numBaseFunctions;
00318       for( ; rit != rend; )
00319       {
00320         F *prit = rit;
00321 
00322         *wit = z * *rit; 
00323         ++rit, ++wit; 
00324 
00325         for( unsigned d = 1; d <= deriv; ++d )
00326         {
00327           #ifndef NDEBUG
00328           const F *const derivEnd = rit + mySize.sizes_[ d ];
00329           #endif
00330           const F *const drend = rit + mySize.sizes_[ d ] - mySize.sizes_[ d-1 ];
00331           for( ; rit != drend ; ++rit, ++wit )
00332             *wit = z * *rit; 
00333           for (unsigned int j=1;j<d;++j) 
00334           {
00335             const F *const drend = rit + mySize.sizes_[ d-j ] - mySize.sizes_[ d-j-1 ];
00336             for( ; rit != drend ; ++prit, ++rit, ++wit )
00337               *wit = F(j) * *prit + z * *rit; 
00338           }
00339           *wit = F(d) * *prit + z * *rit; 
00340           ++prit, ++rit, ++wit;
00341           assert(derivEnd == rit);
00342           rit += size.sizes_[d] - mySize.sizes_[d];
00343           prit += size.sizes_[d-1] - mySize.sizes_[d-1];
00344           const F *const emptyWitEnd = wit + size.sizes_[d] - mySize.sizes_[d];
00345           for ( ; wit != emptyWitEnd; ++wit )
00346             *wit = Zero<F>();
00347         }
00348       }
00349     }
00350   };
00351 
00352 
00353 
00354   // MonomialBasisImpl
00355   // -----------------
00356 
00357   template< class Topology, class F >
00358   class MonomialBasisImpl;
00359 
00360   template< class F >
00361   class MonomialBasisImpl< GenericGeometry::Point, F >
00362   {
00363     typedef MonomialBasisImpl< GenericGeometry::Point, F > This;
00364 
00365   public:
00366     typedef GenericGeometry::Point Topology;
00367     typedef F Field;
00368 
00369     static const unsigned int dimDomain = Topology::dimension;
00370 
00371     typedef FieldVector< Field, dimDomain > DomainVector;
00372 
00373   private:
00374     friend class MonomialBasis< Topology, Field >;
00375     friend class MonomialBasisImpl< GenericGeometry::Prism< Topology >, Field >;
00376     friend class MonomialBasisImpl< GenericGeometry::Pyramid< Topology >, Field >;
00377 
00378     template< int dimD >
00379     void evaluate ( const unsigned int deriv, const unsigned int order,
00380                     const FieldVector< Field, dimD > &x,
00381                     const unsigned int block, const unsigned int *const offsets,
00382                     Field *const values ) const  
00383     {
00384       *values = Unity< F >();
00385       F *const end = values + block;
00386       for( Field *it = values+1 ; it != end; ++it )
00387         *it = Zero< F >();
00388     }
00389 
00390     void integrate ( const unsigned int order,
00391                      const unsigned int *const offsets,
00392                      Field *const values ) const
00393     {
00394       values[ 0 ] = Unity< Field >();
00395     }
00396   };
00397 
00398   template< class BaseTopology, class F >
00399   class MonomialBasisImpl< GenericGeometry::Prism< BaseTopology >, F >
00400   {
00401     typedef MonomialBasisImpl< GenericGeometry::Prism< BaseTopology >, F > This;
00402 
00403   public:
00404     typedef GenericGeometry::Prism< BaseTopology > Topology;
00405     typedef F Field;
00406 
00407     static const unsigned int dimDomain = Topology::dimension;
00408 
00409     typedef FieldVector< Field, dimDomain > DomainVector;
00410 
00411   private:
00412     friend class MonomialBasis< Topology, Field >;
00413     friend class MonomialBasisImpl< GenericGeometry::Prism< Topology >, Field >;
00414     friend class MonomialBasisImpl< GenericGeometry::Pyramid< Topology >, Field >;
00415 
00416     typedef MonomialBasisSize< BaseTopology > BaseSize; 
00417     typedef MonomialBasisSize< Topology > Size; 
00418 
00419     MonomialBasisImpl< BaseTopology, Field > baseBasis_;
00420 
00421     MonomialBasisImpl ()
00422     {}
00423 
00424     template< int dimD >
00425     void evaluate ( const unsigned int deriv, const unsigned int order,
00426                     const FieldVector< Field, dimD > &x,
00427                     const unsigned int block, const unsigned int *const offsets,
00428                     Field *const values ) const
00429     {
00430       typedef MonomialBasisHelper< dimDomain, dimD, Field > Helper;
00431       const BaseSize &size = BaseSize::instance();
00432 
00433       const Field &z = x[ dimDomain-1 ];
00434       
00435       // fill first column
00436       baseBasis_.evaluate( deriv, order, x, block, offsets, values );
00437 
00438       Field *row0 = values;
00439       for( unsigned int k = 1; k <= order; ++k )
00440       {
00441         Field *row1 = values + block*offsets[ k-1 ];
00442         Field *wit = row1 + block*size.sizes_[ k ];
00443         Helper::copy( deriv, wit, row1, k*size.sizes_[ k ], z );
00444         Helper::copy( deriv, wit, row0, size( k-1 ), z );
00445         row0 = row1;
00446       }
00447     }
00448 
00449     void integrate ( const unsigned int order,
00450                      const unsigned int *const offsets,
00451                      Field *const values ) const
00452     {
00453       const BaseSize &size = BaseSize::instance();
00454       const Size &mySize = Size::instance();
00455       // fill first column
00456       baseBasis_.integrate( order, offsets, values );
00457       const unsigned int *const baseSizes = size.sizes_;
00458 
00459       Field *row0 = values;
00460       for( unsigned int k = 1; k <= order; ++k )
00461       {
00462         Field *const row1begin = values + offsets[ k-1 ];
00463         Field *const row1End = row1begin + mySize.sizes_[ k ];
00464         assert( (unsigned int)(row1End - values) <= offsets[ k ] );
00465 
00466         Field *row1 = row1begin;
00467         Field *it = row1begin + baseSizes[ k ];
00468         for( unsigned int j = 1; j <= k; ++j )
00469         {
00470           Field *const end = it + baseSizes[ k ];
00471           assert( (unsigned int)(end - values) <= offsets[ k ] );
00472           for( ; it != end; ++row1, ++it )
00473             *it = (Field( j ) / Field( j+1 )) * (*row1);
00474         }
00475         for( ; it != row1End; ++row0, ++it )
00476           *it = (Field( k ) / Field( k+1 )) * (*row0);
00477         row0 = row1;
00478       }
00479     }
00480   };
00481 
00482   template< class BaseTopology, class F >
00483   class MonomialBasisImpl< GenericGeometry::Pyramid< BaseTopology >, F >
00484   {
00485     typedef MonomialBasisImpl< GenericGeometry::Pyramid< BaseTopology >, F > This;
00486 
00487   public:
00488     typedef GenericGeometry::Pyramid< BaseTopology > Topology;
00489     typedef F Field;
00490 
00491     static const unsigned int dimDomain = Topology::dimension;
00492 
00493     typedef FieldVector< Field, dimDomain > DomainVector;
00494 
00495   private:
00496     friend class MonomialBasis< Topology, Field >;
00497     friend class MonomialBasisImpl< GenericGeometry::Prism< Topology >, Field >;
00498     friend class MonomialBasisImpl< GenericGeometry::Pyramid< Topology >, Field >;
00499 
00500     typedef MonomialBasisSize< BaseTopology > BaseSize; 
00501     typedef MonomialBasisSize< Topology > Size; 
00502 
00503     MonomialBasisImpl< BaseTopology, Field > baseBasis_;
00504 
00505     MonomialBasisImpl ()
00506     {
00507     }
00508 
00509     template< int dimD >
00510     void evaluateSimplexBase ( const unsigned int deriv, const unsigned int order,
00511                                const FieldVector< Field, dimD > &x,
00512                                const unsigned int block, const unsigned int *const offsets,
00513                                Field *const values, 
00514                                const BaseSize &size ) const
00515     {
00516       baseBasis_.evaluate( deriv, order, x, block, offsets, values );
00517     }
00518 
00519     template< int dimD >
00520     void evaluatePyramidBase ( const unsigned int deriv, const unsigned int order,
00521                                const FieldVector< Field, dimD > &x,
00522                                const unsigned int block, const unsigned int *const offsets,
00523                                Field *const values,
00524                                const BaseSize &size ) const
00525     {
00526       Field omz = Unity< Field >() - x[ dimDomain-1 ];
00527 
00528       if( Zero< Field >() < omz ) 
00529       {
00530         const Field invomz = Unity< Field >() / omz;
00531         FieldVector< Field, dimD > y;
00532         for( unsigned int i = 0; i < dimDomain-1; ++i )
00533           y[ i ] = x[ i ] * invomz;
00534       
00535         // fill first column
00536         baseBasis_.evaluate( deriv, order, y, block, offsets, values );
00537 
00538         Field omzk = omz;
00539         for( unsigned int k = 1; k <= order; ++k )
00540         {
00541           Field *it = values + block*offsets[ k-1 ];
00542           Field *const end = it + block*size.sizes_[ k ];
00543           for( ; it != end; ++it )
00544             *it *= omzk;
00545           omzk *= omz;
00546         }
00547       }
00548       else
00549       {
00550         assert( deriv==0 );
00551         *values = Unity< Field >();
00552         for( unsigned int k = 1; k <= order; ++k )
00553         {
00554           Field *it = values + block*offsets[ k-1 ];
00555           Field *const end = it + block*size.sizes_[ k ];
00556           for( ; it != end; ++it )
00557             *it = Zero< Field >();
00558         }
00559       }
00560     }
00561 
00562     template< int dimD >
00563     void evaluate ( const unsigned int deriv, const unsigned int order,
00564                     const FieldVector< Field, dimD > &x,
00565                     const unsigned int block, const unsigned int *const offsets,
00566                     Field *const values ) const
00567     {
00568       typedef MonomialBasisHelper< dimDomain, dimD, Field > Helper;
00569       const BaseSize &size = BaseSize::instance();
00570       
00571       if( GenericGeometry::IsSimplex< Topology >::value )
00572         evaluateSimplexBase( deriv, order, x, block, offsets, values, size );
00573       else
00574         evaluatePyramidBase( deriv, order, x, block, offsets, values, size );
00575 
00576       Field *row0 = values;
00577       for( unsigned int k = 1; k <= order; ++k )
00578       {
00579         Field *row1 = values + block*offsets[ k-1 ];
00580         Field *wit = row1 + block*size.sizes_[ k ];
00581         Helper::copy( deriv, wit, row0, size( k-1 ), x[ dimDomain-1 ] );
00582         row0 = row1;
00583       }
00584     }
00585 
00586     void integrate ( const unsigned int order,
00587                      const unsigned int *const offsets,
00588                      Field *const values ) const
00589     {
00590       const BaseSize &size = BaseSize::instance();
00591 
00592       // fill first column
00593       baseBasis_.integrate( order, offsets, values );
00594 
00595       const unsigned int *const baseSizes = size.sizes_;
00596 
00597       Field *const col0End = values + baseSizes[ 0 ];
00598       for( Field *it = values; it != col0End; ++it )
00599         *it *= Field( 1 ) /  Field( int(dimDomain) );   
00600       Field *row0 = values;
00601 
00602       for( unsigned int k = 1; k <= order; ++k )
00603       {
00604         const Field factor = (Field( 1 ) / Field( k + dimDomain ));
00605 
00606         Field *const row1 = values+offsets[ k-1 ];
00607         Field *const col0End = row1 + baseSizes[ k ];
00608         Field *it = row1;
00609         for( ; it != col0End; ++it )
00610           *it *= factor;
00611         for( unsigned int i = 1; i <= k; ++i )
00612         {
00613           Field *const end = it + baseSizes[ k-i ];
00614           assert( (unsigned int)(end - values) <= offsets[ k ] );
00615           for( ; it != end; ++row0, ++it )
00616             *it = (*row0) * (Field( i ) * factor);
00617         }
00618         row0 = row1;
00619       }
00620     }
00621   };
00622 
00623 
00624 
00625   // MonomialBasis
00626   // -------------
00627 
00628   template< class Topology, class F >
00629   class MonomialBasis
00630   : public MonomialBasisImpl< Topology, F >
00631   {
00632     typedef MonomialBasis< Topology, F > This;
00633     typedef MonomialBasisImpl< Topology, F > Base;
00634 
00635   public:
00636     static const unsigned int dimension = Base::dimDomain;
00637     static const unsigned int dimRange = 1;
00638 
00639     typedef typename Base::Field Field;
00640 
00641     typedef typename Base::DomainVector DomainVector;
00642 
00643     typedef Dune::FieldVector<Field,dimRange> RangeVector;
00644 
00645     typedef MonomialBasisSize<Topology> Size;
00646 
00647     MonomialBasis (unsigned int order)
00648     : Base(),
00649       order_(order),
00650       size_(Size::instance())
00651     {
00652       assert(order<=1024); // avoid wrapping of unsigned int (0-1) order=1024 is quite hight...)
00653     }
00654 
00655     const unsigned int *sizes ( unsigned int order ) const
00656     {
00657       size_.computeSizes( order );
00658       return size_.numBaseFunctions_;
00659     }
00660 
00661     const unsigned int *sizes () const
00662     {
00663       return sizes( order_ );
00664     }
00665 
00666     const unsigned int size () const
00667     {
00668       size_.computeSizes( order_ );
00669       return size_( order_ );
00670     }
00671 
00672     const unsigned int derivSize ( const unsigned int deriv ) const
00673     {
00674       typedef typename GenericGeometry::SimplexTopology< dimension >::type SimplexTopology;
00675       MonomialBasisSize< SimplexTopology >::instance().computeSizes( deriv );
00676       return MonomialBasisSize< SimplexTopology >::instance()( deriv );
00677     }
00678 
00679     const unsigned int order () const
00680     {
00681       return order_ ;
00682     }
00683 
00684     const unsigned int topologyId ( ) const
00685     {
00686       return Topology::id;
00687     }
00688 
00689     void evaluate ( const unsigned int deriv, const DomainVector &x,
00690                     Field *const values ) const
00691     {
00692       Base::evaluate( deriv, order_, x, derivSize( deriv ), sizes( order_ ), values );
00693     }
00694 
00695     template <unsigned int deriv>
00696     void evaluate ( const DomainVector &x,
00697                     Field *const values ) const
00698     {
00699       evaluate( deriv, x, values );
00700     }
00701 
00702     template<unsigned int deriv, class Vector >
00703     void evaluate ( const DomainVector &x,
00704                     Vector &values ) const
00705     {
00706       evaluate<deriv>(x,&(values[0]));
00707     }
00708     template<unsigned int deriv, DerivativeLayout layout >
00709     void evaluate ( const DomainVector &x,
00710                     Derivatives<Field,dimension,1,deriv,layout> *values ) const
00711     {
00712       evaluate<deriv>(x,&(values->block()));
00713     }
00714     template< unsigned int deriv >
00715     void evaluate ( const DomainVector &x,
00716                     FieldVector<Field,Derivatives<Field,dimension,1,deriv,value>::size> *values ) const
00717     {
00718       evaluate(0,x,&(values[0][0]));
00719     }
00720 
00721     template<class Vector >
00722     void evaluate ( const DomainVector &x,
00723                     Vector &values ) const
00724     {
00725       evaluate<0>(x,&(values[0]));
00726     }
00727 
00728     template< class DVector, class RVector >
00729     void evaluate ( const DVector &x, RVector &values ) const
00730     {
00731       assert( DVector::size == dimension);
00732       DomainVector bx;
00733       for( int d = 0; d < dimension; ++d )
00734         field_cast( x[ d ], bx[ d ] );
00735       evaluate<0>( bx, values );
00736     }
00737 
00738     void integrate ( Field *const values ) const
00739     {
00740       Base::integrate( order_, sizes( order_ ), values );
00741     }
00742     template <class Vector>
00743     void integrate ( Vector &values ) const
00744     {
00745       integrate( &(values[ 0 ]) );
00746     }
00747   private:
00748     MonomialBasis(const This&);
00749     This& operator=(const This&);
00750     unsigned int order_;
00751     Size &size_; 
00752   };
00753 
00754 
00755 
00756   // StdMonomialBasis
00757   // ----------------
00758 
00759   template< int dim,class F >
00760   class StandardMonomialBasis
00761   : public MonomialBasis< typename GenericGeometry::SimplexTopology< dim >::type, F >
00762   {
00763     typedef StandardMonomialBasis< dim, F > This;
00764     typedef MonomialBasis< typename GenericGeometry::SimplexTopology< dim >::type, F > Base;
00765 
00766   public:
00767     typedef typename GenericGeometry::SimplexTopology< dim >::type Topology;
00768     static const int dimension = dim;
00769 
00770     StandardMonomialBasis ( unsigned int order )
00771     : Base( order )
00772     {}
00773   };
00774 
00775 
00776 
00777   // StandardBiMonomialBasis
00778   // -----------------------
00779 
00780   template< int dim, class F >
00781   class StandardBiMonomialBasis
00782   : public MonomialBasis< typename GenericGeometry::CubeTopology< dim >::type, F >
00783   {
00784     typedef StandardBiMonomialBasis< dim, F > This;
00785     typedef MonomialBasis< typename GenericGeometry::CubeTopology< dim >::type, F > Base;
00786 
00787   public:
00788     typedef typename GenericGeometry::CubeTopology< dim >::type Topology;
00789     static const int dimension = dim;
00790 
00791     StandardBiMonomialBasis ( unsigned int order )
00792     : Base( order )
00793     {}
00794   };
00795 
00796   // -----------------------------------------------------------
00797   // -----------------------------------------------------------
00798   // VirtualMonomialBasis
00799   // -------------------
00800 
00801   template< int dim, class F >
00802   class VirtualMonomialBasis 
00803   {
00804     typedef VirtualMonomialBasis< dim, F > This;
00805 
00806   public:
00807     typedef F Field;
00808     typedef F StorageField;
00809     static const int dimension = dim;
00810     static const unsigned int dimRange = 1;
00811 
00812     typedef FieldVector<Field,dimension> DomainVector;
00813     typedef FieldVector<Field,dimRange> RangeVector;
00814 
00815     explicit VirtualMonomialBasis(unsigned int topologyId,
00816                                   unsigned int order)
00817       : order_(order), topologyId_(topologyId) {}
00818 
00819     virtual ~VirtualMonomialBasis() {}
00820 
00821     virtual const unsigned int *sizes ( ) const = 0;
00822 
00823     const unsigned int size ( ) const
00824     {
00825       return sizes( )[ order_ ];
00826     }
00827 
00828     const unsigned int order () const
00829     {
00830       return order_;
00831     }
00832 
00833     const unsigned int topologyId ( ) const
00834     {
00835       return topologyId_;
00836     }
00837 
00838     virtual void evaluate ( const unsigned int deriv, const DomainVector &x,
00839                             Field *const values ) const = 0;
00840     template < unsigned int deriv >
00841     void evaluate ( const DomainVector &x,
00842                     Field *const values ) const 
00843     {
00844       evaluate( deriv, x, values );
00845     }
00846     template < unsigned int deriv, int size >
00847     void evaluate ( const DomainVector &x,
00848                     Dune::FieldVector<Field,size> *const values ) const
00849     {
00850       evaluate( deriv, x, &(values[0][0]) );
00851     }
00852     template<unsigned int deriv, DerivativeLayout layout >
00853     void evaluate ( const DomainVector &x,
00854                     Derivatives<Field,dimension,1,deriv,layout> *values ) const
00855     {
00856       evaluate<deriv>(x,&(values->block()));
00857     }
00858     template <unsigned int deriv, class Vector>
00859     void evaluate ( const DomainVector &x,
00860                     Vector &values ) const
00861     {
00862       evaluate<deriv>( x, &(values[ 0 ]) );
00863     }
00864     template< class Vector >
00865     void evaluate ( const DomainVector &x,
00866                     Vector &values ) const
00867     {
00868       evaluate<0>(x,values);
00869     }
00870     template< class DVector, class RVector >
00871     void evaluate ( const DVector &x, RVector &values ) const
00872     {
00873       assert( DVector::size == dimension);
00874       DomainVector bx;
00875       for( int d = 0; d < dimension; ++d )
00876         field_cast( x[ d ], bx[ d ] );
00877       evaluate<0>( bx, values );
00878     }
00879     template< unsigned int deriv, class DVector, class RVector >
00880     void evaluate ( const DVector &x, RVector &values ) const
00881     {
00882       assert( DVector::size == dimension);
00883       DomainVector bx;
00884       for( int d = 0; d < dimension; ++d )
00885         field_cast( x[ d ], bx[ d ] );
00886       evaluate<deriv>( bx, values );
00887     }
00888 
00889     virtual void integrate ( Field *const values ) const = 0;
00890     template <class Vector>
00891     void integrate ( Vector &values ) const
00892     {
00893       integrate( &(values[ 0 ]) );
00894     }
00895   protected:
00896     unsigned int order_;
00897     unsigned int topologyId_;
00898   };
00899 
00900   template< class Topology, class F >
00901   class VirtualMonomialBasisImpl 
00902   : public VirtualMonomialBasis< Topology::dimension, F >
00903   {
00904     typedef VirtualMonomialBasis< Topology::dimension, F > Base;
00905     typedef VirtualMonomialBasisImpl< Topology, F > This;
00906 
00907   public:
00908     typedef typename Base::Field Field;
00909     typedef typename Base::DomainVector DomainVector;
00910 
00911     VirtualMonomialBasisImpl(unsigned int order)
00912     : Base(Topology::id,order), basis_(order) 
00913     {}
00914 
00915     const unsigned int *sizes ( ) const
00916     {
00917       return basis_.sizes(order_);
00918     }
00919 
00920     void evaluate ( const unsigned int deriv, const DomainVector &x,
00921                     Field *const values ) const
00922     {
00923       basis_.evaluate(deriv,x,values);
00924     }
00925 
00926     void integrate ( Field *const values ) const
00927     {
00928       basis_.integrate(values);
00929     }
00930 
00931   private:
00932     MonomialBasis<Topology,Field> basis_;
00933     using Base::order_;
00934   };
00935 
00936   // MonomialBasisFactory
00937   // --------------------
00938   
00939   template< int dim, class F >
00940   struct MonomialBasisFactory;
00941 
00942   template< int dim, class F >
00943   struct MonomialBasisFactoryTraits
00944   {
00945     static const unsigned int dimension = dim;
00946     typedef unsigned int Key;
00947     typedef const VirtualMonomialBasis< dimension, F > Object;
00948     typedef MonomialBasisFactory<dim,F> Factory;
00949   };
00950 
00951   template< int dim, class F >
00952   struct MonomialBasisFactory : 
00953     public TopologyFactory< MonomialBasisFactoryTraits<dim,F> >
00954   {
00955     static const unsigned int dimension = dim;
00956     typedef F StorageField;
00957     typedef MonomialBasisFactoryTraits<dim,F> Traits;
00958 
00959     typedef typename Traits::Key Key;
00960     typedef typename Traits::Object Object;
00961 
00962     template < int dd, class FF >
00963     struct EvaluationBasisFactory
00964     {
00965       typedef MonomialBasisFactory<dd,FF> Type;
00966     };
00967 
00968     template< class Topology >
00969     static Object* createObject ( const Key &order )
00970     {
00971       return new VirtualMonomialBasisImpl< Topology, StorageField >( order );
00972     }
00973   };
00974 
00975 
00976 
00977   // MonomialBasisProvider
00978   // ---------------------
00979 
00980   template< int dim, class SF >
00981   struct MonomialBasisProvider 
00982   : public TopologySingletonFactory< MonomialBasisFactory< dim, SF > >
00983   {
00984     static const unsigned int dimension = dim;
00985     typedef SF StorageField;
00986     template < int dd, class FF >
00987     struct EvaluationBasisFactory
00988     {
00989       typedef MonomialBasisProvider<dd,FF> Type;
00990     };
00991   };
00992 
00993 }
00994 
00995 #endif