dune-pdelab  2.4-dev
localbasiscache.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil -*-
2 #ifndef DUNE_PDELAB_LOCALBASISCACHE_HH
3 #define DUNE_PDELAB_LOCALBASISCACHE_HH
4 
5 #include<vector>
6 #include<map>
7 
8 #include<dune/common/exceptions.hh>
9 
10 namespace Dune {
11  namespace PDELab {
12 
14  template<class LocalBasisType>
16  {
17  typedef typename LocalBasisType::Traits::DomainFieldType DomainFieldType;
18  typedef typename LocalBasisType::Traits::DomainType DomainType;
19  typedef typename LocalBasisType::Traits::RangeType RangeType;
20  typedef typename LocalBasisType::Traits::JacobianType JacobianType;
21 
22  struct less_than
23  {
24  bool operator() (const DomainType& v1, const DomainType& v2) const
25  {
26  for (typename DomainType::size_type i=0; i<DomainType::dimension; i++)
27  {
28  if ( v1[i] < v2[i]-1e-5 ) return true; // is less than
29  if ( v1[i] > v2[i]+1e-5 ) return false; // is greater than
30  }
31  return false; // is equal
32  }
33  };
34 
35  typedef std::map<DomainType,std::vector<RangeType>,less_than> FunctionCache;
36  typedef std::map<DomainType,std::vector<JacobianType>,less_than> JacobianCache;
37 
38  public:
39 
42 
44  const std::vector<RangeType>&
45  evaluateFunction (const DomainType& position, const LocalBasisType& localbasis) const
46  {
47  typename FunctionCache::iterator it = functioncache.find(position);
48  if (it!=functioncache.end()) return it->second;
49  std::vector<RangeType> values;
50  localbasis.evaluateFunction(position,values);
51  it = functioncache.insert(functioncache.begin(),std::pair<DomainType,std::vector<RangeType> >(position,values));
52  return it->second;
53  }
54 
56  const std::vector<JacobianType>&
57  evaluateJacobian (const DomainType& position, const LocalBasisType& localbasis) const
58  {
59  typename JacobianCache::iterator it = jacobiancache.find(position);
60  if (it!=jacobiancache.end()) return it->second;
61  std::vector<JacobianType> values;
62  localbasis.evaluateJacobian(position,values);
63  it = jacobiancache.insert(jacobiancache.begin(),std::pair<DomainType,std::vector<JacobianType> >(position,values));
64  return it->second;
65  }
66 
67  private:
68  mutable FunctionCache functioncache;
69  mutable JacobianCache jacobiancache;
70  };
71 
72  }
73 }
74 
75 #endif
LocalBasisCache()
constructor
Definition: localbasiscache.hh:41
const E & e
Definition: interpolate.hh:172
Definition: adaptivity.hh:27
store values of basis functions and gradients in a cache
Definition: localbasiscache.hh:15
const std::vector< RangeType > & evaluateFunction(const DomainType &position, const LocalBasisType &localbasis) const
evaluate basis functions at a point
Definition: localbasiscache.hh:45
const std::vector< JacobianType > & evaluateJacobian(const DomainType &position, const LocalBasisType &localbasis) const
evaluate Jacobians at a point
Definition: localbasiscache.hh:57