00001 #ifndef DUNE_RAVIARTTHOMASINTERPOLATION_HH
00002 #define DUNE_RAVIARTTHOMASINTERPOLATION_HH
00003
00004 #include <fstream>
00005 #include <utility>
00006
00007 #include <dune/common/exceptions.hh>
00008 #include <dune/common/forloop.hh>
00009
00010 #include <dune/grid/common/quadraturerules/gaussquadrature.hh>
00011 #include <dune/grid/genericgeometry/referencemappings.hh>
00012 #include <dune/grid/genericgeometry/referenceelements.hh>
00013
00014 #include <dune/localfunctions/common/localkey.hh>
00015 #include <dune/grid/common/topologyfactory.hh>
00016 #include <dune/localfunctions/utility/interpolationhelper.hh>
00017 #include <dune/localfunctions/utility/lfematrix.hh>
00018 #include <dune/localfunctions/utility/polynomialbasis.hh>
00019 #include <dune/localfunctions/orthonormal/orthonormalbasis.hh>
00020
00021 namespace Dune
00022 {
00023
00024
00025
00026 template < unsigned int dim >
00027 struct RaviartThomasCoefficientsFactory;
00028 template < unsigned int dim, class Field >
00029 struct RaviartThomasL2InterpolationFactory;
00030
00031 class LocalCoefficientsContainer
00032 {
00033 typedef LocalCoefficientsContainer This;
00034
00035 public:
00036 template <class Setter>
00037 LocalCoefficientsContainer ( const Setter &setter )
00038 {
00039 setter.setLocalKeys(localKey_);
00040 }
00041
00042 const LocalKey &localKey ( const unsigned int i ) const
00043 {
00044 assert( i < size() );
00045 return localKey_[ i ];
00046 }
00047
00048 unsigned int size () const
00049 {
00050 return localKey_.size();
00051 }
00052
00053 private:
00054 std::vector< LocalKey > localKey_;
00055 };
00056
00057 template < unsigned int dim >
00058 struct RaviartThomasCoefficientsFactoryTraits
00059 {
00060 static const unsigned int dimension = dim;
00061 typedef const LocalCoefficientsContainer Object;
00062 typedef unsigned int Key;
00063 typedef RaviartThomasCoefficientsFactory<dim> Factory;
00064 };
00065 template < unsigned int dim >
00066 struct RaviartThomasCoefficientsFactory :
00067 public TopologyFactory< RaviartThomasCoefficientsFactoryTraits<dim> >
00068 {
00069 typedef RaviartThomasCoefficientsFactoryTraits<dim> Traits;
00070 template <class Topology>
00071 static typename Traits::Object *createObject( const typename Traits::Key &key )
00072 {
00073 typedef RaviartThomasL2InterpolationFactory<dim,double> InterpolationFactory;
00074 if (! supports<Topology>(key) )
00075 return 0;
00076 typename InterpolationFactory::Object *interpol
00077 = InterpolationFactory::template create<Topology>(key);
00078 typename Traits::Object *localKeys = new typename Traits::Object(*interpol);
00079 InterpolationFactory::release(interpol);
00080 return localKeys;
00081 }
00082 template< class Topology >
00083 static bool supports ( const typename Traits::Key &key )
00084 {
00085 return GenericGeometry::IsSimplex<Topology>::value;
00086 }
00087 };
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101 template <unsigned int dim, class Field>
00102 struct RTL2InterpolationBuilder
00103 {
00104 static const unsigned int dimension = dim;
00105 typedef OrthonormalBasisFactory<dimension,Field> TestBasisFactory;
00106 typedef typename TestBasisFactory::Object TestBasis;
00107 typedef OrthonormalBasisFactory<dimension-1,Field> TestFaceBasisFactory;
00108 typedef typename TestFaceBasisFactory::Object TestFaceBasis;
00109 typedef FieldVector<Field,dimension> Normal;
00110
00111 RTL2InterpolationBuilder()
00112 {}
00113
00114 ~RTL2InterpolationBuilder()
00115 {
00116 TestBasisFactory::release(testBasis_);
00117 for (unsigned int i=0;i<faceStructure_.size();++i)
00118 TestFaceBasisFactory::release(faceStructure_[i].basis_);
00119 }
00120
00121 unsigned int topologyId() const
00122 {
00123 return topologyId_;
00124 }
00125 unsigned int order() const
00126 {
00127 return order_;
00128 }
00129 unsigned int faceSize() const
00130 {
00131 return faceSize_;
00132 }
00133 TestBasis *testBasis() const
00134 {
00135 return testBasis_;
00136 }
00137 TestFaceBasis *testFaceBasis( unsigned int f ) const
00138 {
00139 assert( f < faceSize() );
00140 return faceStructure_[f].basis_;
00141 }
00142 const Normal &normal( unsigned int f ) const
00143 {
00144 return *(faceStructure_[f].normal_);
00145 }
00146
00147 template <class Topology>
00148 void build(unsigned int order)
00149 {
00150 order_ = order;
00151 topologyId_ = Topology::id;
00152 if (order>0)
00153 testBasis_ = TestBasisFactory::template create<Topology>(order-1);
00154 else
00155 testBasis_ = 0;
00156 const unsigned int size = GenericGeometry::Size<Topology,1>::value;
00157 faceSize_ = size;
00158 faceStructure_.reserve( faceSize_ );
00159 ForLoop< Creator<Topology>::template GetCodim,0,size-1>::apply(order, faceStructure_ );
00160 assert( faceStructure_.size() == faceSize_ );
00161 }
00162
00163 private:
00164 struct FaceStructure
00165 {
00166 FaceStructure( TestFaceBasis *tfb,
00167 const Normal &n )
00168 : basis_(tfb), normal_(&n)
00169 {
00170 }
00171 TestFaceBasis *basis_;
00172 const Dune::FieldVector<Field,dimension> *normal_;
00173 };
00174 template < class Topology >
00175 struct Creator
00176 {
00177 template < int face >
00178 struct GetCodim
00179 {
00180 typedef typename GenericGeometry::SubTopology<Topology,1,face>::type FaceTopology;
00181 static void apply( const unsigned int order,
00182 std::vector<FaceStructure> &faceStructure )
00183 {
00184 faceStructure.push_back(
00185 FaceStructure(
00186 TestFaceBasisFactory::template create<FaceTopology>(order),
00187 GenericGeometry::ReferenceElement<Topology,Field>::integrationOuterNormal(face)
00188 ) );
00189 }
00190 };
00191 };
00192
00193 std::vector<FaceStructure> faceStructure_;
00194 TestBasis *testBasis_;
00195 unsigned int topologyId_, order_, faceSize_;
00196 };
00197
00198
00199
00200 template< unsigned int dimension, class F>
00201 class RaviartThomasL2Interpolation
00202 : public InterpolationHelper<F,dimension>
00203 {
00204 typedef RaviartThomasL2Interpolation< dimension, F > This;
00205 typedef InterpolationHelper<F,dimension> Base;
00206
00207 public:
00208 typedef F Field;
00209 typedef RTL2InterpolationBuilder<dimension,Field> Builder;
00210 RaviartThomasL2Interpolation( )
00211 : order_(0),
00212 size_(0)
00213 {
00214 }
00215
00216 template< class Function, class Fy >
00217 void interpolate ( const Function &function, std::vector< Fy > &coefficients ) const
00218 {
00219 coefficients.resize(size());
00220 typename Base::template Helper<Function,std::vector<Fy>,true> func( function,coefficients );
00221 interpolate(func);
00222 }
00223 template< class Basis, class Matrix >
00224 void interpolate ( const Basis &basis, Matrix &matrix ) const
00225 {
00226 matrix.resize( size(), basis.size() );
00227 typename Base::template Helper<Basis,Matrix,false> func( basis,matrix );
00228 interpolate(func);
00229 }
00230
00231 unsigned int order() const
00232 {
00233 return order_;
00234 }
00235 unsigned int size() const
00236 {
00237 return size_;
00238 }
00239 template <class Topology>
00240 void build( unsigned int order )
00241 {
00242 size_ = 0;
00243 order_ = order;
00244 builder_.template build<Topology>(order_);
00245 if (builder_.testBasis())
00246 size_ += dimension*builder_.testBasis()->size();
00247 for ( unsigned int f=0;f<builder_.faceSize();++f )
00248 if (builder_.testFaceBasis(f))
00249 size_ += builder_.testFaceBasis(f)->size();
00250 }
00251
00252 void setLocalKeys(std::vector< LocalKey > &keys) const
00253 {
00254 keys.resize(size());
00255 unsigned int row = 0;
00256 for (unsigned int f=0;f<builder_.faceSize();++f)
00257 {
00258 if (builder_.faceSize())
00259 for (unsigned int i=0;i<builder_.testFaceBasis(f)->size();++i,++row)
00260 keys[row] = LocalKey(f,1,i);
00261 }
00262 if (builder_.testBasis())
00263 for (unsigned int i=0;i<builder_.testBasis()->size()*dimension;++i,++row)
00264 keys[row] = LocalKey(0,0,i);
00265 assert( row == size() );
00266 }
00267
00268 protected:
00269 template< class Func, class Container, bool type >
00270 void interpolate ( typename Base::template Helper<Func,Container,type> &func ) const
00271 {
00272 const unsigned int topologyId = builder_.topologyId();
00273
00274 std::vector< Field > testBasisVal;
00275
00276 for (unsigned int i=0;i<size();++i)
00277 for (unsigned int j=0;j<func.size();++j)
00278 func.set(i,j,0);
00279
00280 unsigned int row = 0;
00281
00282
00283 typedef Dune::GenericGeometry::GaussQuadratureProvider< dimension-1, Field >
00284 FaceQuadratureProvider;
00285 typedef Dune::GenericGeometry::ReferenceMappings< Field, dimension > RefMappings;
00286 typedef typename RefMappings::Container::
00287 template Codim< 1 >::Mapping Mapping;
00288
00289 for (unsigned int f=0;f<builder_.faceSize();++f)
00290 {
00291 if (!builder_.testFaceBasis(f))
00292 continue;
00293 testBasisVal.resize(builder_.testFaceBasis(f)->size());
00294
00295 const Mapping& mapping = RefMappings::container( topologyId ).
00296 template mapping< 1 >( f );
00297 const unsigned int subTopologyId = mapping.topologyId();
00298 const typename FaceQuadratureProvider::Object *faceQuad = FaceQuadratureProvider::create( subTopologyId, 2*order_+2 );
00299
00300 const unsigned int quadratureSize = faceQuad->size();
00301 for( unsigned int qi = 0; qi < quadratureSize; ++qi )
00302 {
00303 builder_.testFaceBasis(f)->template evaluate<0>(faceQuad->position(qi),testBasisVal);
00304 fillBnd( row, testBasisVal,
00305 func.evaluate( mapping.global( faceQuad->position(qi) ) ),
00306 builder_.normal(f), faceQuad->weight(qi),
00307 func);
00308 }
00309
00310 FaceQuadratureProvider::release( faceQuad );
00311
00312 row += builder_.testFaceBasis(f)->size();
00313 }
00314
00315 if (builder_.testBasis())
00316 {
00317 testBasisVal.resize(builder_.testBasis()->size());
00318
00319 typedef Dune::GenericGeometry::GaussQuadratureProvider< dimension, Field > QuadratureProvider;
00320 const typename QuadratureProvider::Object *elemQuad = QuadratureProvider::create(topologyId,2*order_+1);
00321
00322 const unsigned int quadratureSize = elemQuad->size();
00323 for( unsigned int qi = 0; qi < quadratureSize; ++qi )
00324 {
00325 builder_.testBasis()->template evaluate<0>(elemQuad->position(qi),testBasisVal);
00326 fillInterior( row, testBasisVal,
00327 func.evaluate(elemQuad->position(qi)),
00328 elemQuad->weight(qi),
00329 func );
00330 }
00331
00332 QuadratureProvider::release( elemQuad );
00333
00334 row += builder_.testBasis()->size()*dimension;
00335 }
00336 assert(row==size());
00337 }
00338
00339 private:
00341 template <class MVal, class RTVal,class Matrix>
00342 void fillBnd (unsigned int startRow,
00343 const MVal &mVal,
00344 const RTVal &rtVal,
00345 const FieldVector<Field,dimension> &normal,
00346 const Field &weight,
00347 Matrix &matrix) const
00348 {
00349 const unsigned int endRow = startRow+mVal.size();
00350 typename RTVal::const_iterator rtiter = rtVal.begin();
00351 for ( unsigned int col = 0; col < rtVal.size() ; ++rtiter,++col)
00352 {
00353 Field cFactor = (*rtiter)*normal;
00354 typename MVal::const_iterator miter = mVal.begin();
00355 for (unsigned int row = startRow;
00356 row!=endRow; ++miter, ++row )
00357 {
00358 matrix.add(row,col, (weight*cFactor)*(*miter) );
00359 }
00360 assert( miter == mVal.end() );
00361 }
00362 }
00363 template <class MVal, class RTVal,class Matrix>
00364 void fillInterior (unsigned int startRow,
00365 const MVal &mVal,
00366 const RTVal &rtVal,
00367 Field weight,
00368 Matrix &matrix) const
00369 {
00370 const unsigned int endRow = startRow+mVal.size()*dimension;
00371 typename RTVal::const_iterator rtiter = rtVal.begin();
00372 for ( unsigned int col = 0; col < rtVal.size() ; ++rtiter,++col)
00373 {
00374 typename MVal::const_iterator miter = mVal.begin();
00375 for (unsigned int row = startRow;
00376 row!=endRow; ++miter,row+=dimension )
00377 {
00378 for (unsigned int i=0;i<dimension;++i)
00379 {
00380 matrix.add(row+i,col, (weight*(*miter))*(*rtiter)[i] );
00381 }
00382 }
00383 assert( miter == mVal.end() );
00384 }
00385 }
00386
00387 Builder builder_;
00388 unsigned int order_;
00389 unsigned int size_;
00390 };
00391
00392 template < unsigned int dim, class F >
00393 struct RaviartThomasL2InterpolationFactoryTraits
00394 {
00395 static const unsigned int dimension = dim;
00396 typedef unsigned int Key;
00397 typedef const RaviartThomasL2Interpolation<dim,F> Object;
00398 typedef RaviartThomasL2InterpolationFactory<dim,F> Factory;
00399 };
00400 template < unsigned int dim, class Field >
00401 struct RaviartThomasL2InterpolationFactory :
00402 public TopologyFactory< RaviartThomasL2InterpolationFactoryTraits<dim,Field> >
00403 {
00404 typedef RaviartThomasL2InterpolationFactoryTraits<dim,Field> Traits;
00405 typedef RTL2InterpolationBuilder<dim,Field> Builder;
00406 typedef typename Traits::Object Object;
00407 typedef typename remove_const<Object>::type NonConstObject;
00408 template <class Topology>
00409 static typename Traits::Object *createObject( const typename Traits::Key &key )
00410 {
00411 if ( !supports<Topology>(key) )
00412 return 0;
00413 NonConstObject *interpol = new NonConstObject();
00414 interpol->template build<Topology>(key);
00415 return interpol;
00416 }
00417 template< class Topology >
00418 static bool supports ( const typename Traits::Key &key )
00419 {
00420 return GenericGeometry::IsSimplex<Topology>::value;
00421 }
00422 };
00423 }
00424 #endif // DUNE_RAVIARTTHOMASINTERPOLATION_HH
00425