Skip to content
Snippets Groups Projects
Commit ad369fcd authored by graeser's avatar graeser Committed by graeser@PCPOOL.MI.FU-BERLIN.DE
Browse files

Use new virtual interface if DUNE_VIRT* is not enabled

[[Imported from SVN: r3131]]
parent 1509bedb
No related branches found
No related tags found
No related merge requests found
...@@ -8,6 +8,7 @@ ...@@ -8,6 +8,7 @@
#include <dune/grid/common/genericreferenceelements.hh> #include <dune/grid/common/genericreferenceelements.hh>
#include <dune/localfunctions/common/virtualinterface.hh>
#include <dune/localfunctions/p1.hh> #include <dune/localfunctions/p1.hh>
#include <dune/localfunctions/q1.hh> #include <dune/localfunctions/q1.hh>
#include <dune/localfunctions/prismp1.hh> #include <dune/localfunctions/prismp1.hh>
...@@ -31,23 +32,50 @@ class GenericMultigridTransfer { ...@@ -31,23 +32,50 @@ class GenericMultigridTransfer {
template <class DT, class RT, int dim> template <class DT, class RT, int dim>
struct P1ElementFactory { struct P1ElementFactory {
const Dune::LocalFiniteElementInterface<DT, RT, dim>* get(const Dune::GeometryType& type) { #ifdef DUNE_VIRTUAL_SHAPEFUNCTIONS
if (type.isSimplex()) typedef Dune::LocalFiniteElementInterface<DT, RT, dim> type;
const type* get(const Dune::GeometryType& gt) {
if (gt.isSimplex())
return &simplexBaseSet_; return &simplexBaseSet_;
else if (type.isCube()) else if (gt.isCube())
return &cubeBaseSet_; return &cubeBaseSet_;
else if (type.isPrism()) else if (gt.isPrism())
// This cast is necessary because otherwise the code wouldn't compile for dim!=3 // This cast is necessary because otherwise the code wouldn't compile for dim!=3
return (Dune::LocalFiniteElementInterface<DT, RT, dim>*)&prismBaseSet_; return (Dune::LocalFiniteElementInterface<DT, RT, dim>*)&prismBaseSet_;
else else
DUNE_THROW(Dune::NotImplemented, "transfer operators for " << type); DUNE_THROW(Dune::NotImplemented, "transfer operators for " << gt);
} }
private: private:
Dune::P1LocalFiniteElement<DT, RT, dim> simplexBaseSet_; Dune::P1LocalFiniteElement<DT, RT, dim> simplexBaseSet_;
Dune::Q1LocalFiniteElement<DT, RT, dim> cubeBaseSet_; Dune::Q1LocalFiniteElement<DT, RT, dim> cubeBaseSet_;
Dune::PrismP1LocalFiniteElement<DT, RT> prismBaseSet_; Dune::PrismP1LocalFiniteElement<DT, RT> prismBaseSet_;
#else
private:
typedef typename Dune::P1LocalFiniteElement<DT, RT, dim>::Traits LocalBasisTraits;
public:
typedef Dune::C0LocalFiniteElementVirtualInterface<LocalBasisTraits> type;
const type* get(const Dune::GeometryType& gt) {
if (gt.isSimplex())
return &simplexBaseSet_;
else if (gt.isCube())
return &cubeBaseSet_;
else if (gt.isPrism())
// This cast is necessary because otherwise the code wouldn't compile for dim!=3
return (type*)&prismBaseSet_;
else
DUNE_THROW(Dune::NotImplemented, "transfer operators for " << gt);
}
private:
Dune::C0LocalFiniteElementVirtualImp<Dune::P1LocalFiniteElement<DT, RT, dim> > simplexBaseSet_;
Dune::C0LocalFiniteElementVirtualImp<Dune::Q1LocalFiniteElement<DT, RT, dim> > cubeBaseSet_;
Dune::C0LocalFiniteElementVirtualImp<Dune::PrismP1LocalFiniteElement<DT, RT> > prismBaseSet_;
#endif
}; };
...@@ -110,6 +138,7 @@ public: ...@@ -110,6 +138,7 @@ public:
// A factory for the shape functions // A factory for the shape functions
P1ElementFactory<ctype,field_type,dim> p1ElementFactory; P1ElementFactory<ctype,field_type,dim> p1ElementFactory;
typedef typename P1ElementFactory<ctype,field_type,dim>::type FEType;
matrix.setSize(rows,cols); matrix.setSize(rows,cols);
matrix.setBuildMode(TransferOperatorType::random); matrix.setBuildMode(TransferOperatorType::random);
...@@ -131,7 +160,8 @@ public: ...@@ -131,7 +160,8 @@ public:
typedef typename EntityType::HierarchicIterator HierarchicIterator; typedef typename EntityType::HierarchicIterator HierarchicIterator;
// Get local finite element // Get local finite element
const Dune::LocalFiniteElementInterface<ctype, field_type, dim>* coarseBaseSet = p1ElementFactory.get(cIt->type()); // const Dune::LocalFiniteElementInterface<ctype, field_type, dim>* coarseBaseSet = p1ElementFactory.get(cIt->type());
const FEType* coarseBaseSet = p1ElementFactory.get(cIt->type());
const int numCoarseBaseFct = coarseBaseSet->localBasis().size(); const int numCoarseBaseFct = coarseBaseSet->localBasis().size();
...@@ -152,7 +182,8 @@ public: ...@@ -152,7 +182,8 @@ public:
const typename EntityType::Geometry& fGeometryInFather = fIt->geometryInFather(); const typename EntityType::Geometry& fGeometryInFather = fIt->geometryInFather();
// Get local finite element // Get local finite element
const Dune::LocalFiniteElementInterface<ctype, field_type, dim>* fineBaseSet = p1ElementFactory.get(fIt->type());; // const Dune::LocalFiniteElementInterface<ctype, field_type, dim>* fineBaseSet = p1ElementFactory.get(fIt->type());;
const FEType* fineBaseSet = p1ElementFactory.get(fIt->type());;
const int numFineBaseFct = fineBaseSet->localBasis().size(); const int numFineBaseFct = fineBaseSet->localBasis().size();
...@@ -193,7 +224,8 @@ public: ...@@ -193,7 +224,8 @@ public:
for (; cIt != cEndIt; ++cIt) { for (; cIt != cEndIt; ++cIt) {
// Get local finite element // Get local finite element
const Dune::LocalFiniteElementInterface<ctype, field_type, dim>* coarseBaseSet = p1ElementFactory.get(cIt->type()); // const Dune::LocalFiniteElementInterface<ctype, field_type, dim>* coarseBaseSet = p1ElementFactory.get(cIt->type());
const FEType* coarseBaseSet = p1ElementFactory.get(cIt->type());
const int numCoarseBaseFct = coarseBaseSet->localBasis().size(); const int numCoarseBaseFct = coarseBaseSet->localBasis().size();
...@@ -217,7 +249,8 @@ public: ...@@ -217,7 +249,8 @@ public:
const typename EntityType::Geometry& fGeometryInFather = fIt->geometryInFather(); const typename EntityType::Geometry& fGeometryInFather = fIt->geometryInFather();
// Get local finite element // Get local finite element
const Dune::LocalFiniteElementInterface<ctype, field_type, dim>* fineBaseSet = p1ElementFactory.get(fIt->type()); // const Dune::LocalFiniteElementInterface<ctype, field_type, dim>* fineBaseSet = p1ElementFactory.get(fIt->type());
const FEType* fineBaseSet = p1ElementFactory.get(fIt->type());
const int numFineBaseFct = fineBaseSet->localBasis().size(); const int numFineBaseFct = fineBaseSet->localBasis().size();
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment