DUNE PDELab (2.7)

localbasiscache.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4#ifndef DUNE_PDELAB_FINITEELEMENT_LOCALBASISCACHE_HH
5#define DUNE_PDELAB_FINITEELEMENT_LOCALBASISCACHE_HH
6
7#include<vector>
8#include<map>
9
11
12namespace Dune {
13 namespace PDELab {
14
16 template<class LocalBasisType>
18 {
19 typedef typename LocalBasisType::Traits::DomainFieldType DomainFieldType;
20 typedef typename LocalBasisType::Traits::DomainType DomainType;
21 typedef typename LocalBasisType::Traits::RangeType RangeType;
22 typedef typename LocalBasisType::Traits::JacobianType JacobianType;
23
24 struct less_than
25 {
26 bool operator() (const DomainType& v1, const DomainType& v2) const
27 {
28 for (typename DomainType::size_type i=0; i<DomainType::dimension; i++)
29 {
30 if ( v1[i] < v2[i]-1e-5 ) return true; // is less than
31 if ( v1[i] > v2[i]+1e-5 ) return false; // is greater than
32 }
33 return false; // is equal
34 }
35 };
36
37 typedef std::map<DomainType,std::vector<RangeType>,less_than> FunctionCache;
38 typedef std::map<DomainType,std::vector<JacobianType>,less_than> JacobianCache;
39
40 public:
41
44
46 const std::vector<RangeType>&
47 evaluateFunction (const DomainType& position, const LocalBasisType& localbasis) const
48 {
49 typename FunctionCache::iterator it = functioncache.find(position);
50 if (it!=functioncache.end()) return it->second;
51 std::vector<RangeType> values;
52 localbasis.evaluateFunction(position,values);
53 it = functioncache.insert(functioncache.begin(),std::pair<DomainType,std::vector<RangeType> >(position,values));
54 return it->second;
55 }
56
58 const std::vector<JacobianType>&
59 evaluateJacobian (const DomainType& position, const LocalBasisType& localbasis) const
60 {
61 typename JacobianCache::iterator it = jacobiancache.find(position);
62 if (it!=jacobiancache.end()) return it->second;
63 std::vector<JacobianType> values;
64 localbasis.evaluateJacobian(position,values);
65 it = jacobiancache.insert(jacobiancache.begin(),std::pair<DomainType,std::vector<JacobianType> >(position,values));
66 return it->second;
67 }
68
69 private:
70 mutable FunctionCache functioncache;
71 mutable JacobianCache jacobiancache;
72 };
73
74 }
75}
76
77#endif // DUNE_PDELAB_FINITEELEMENT_LOCALBASISCACHE_HH
store values of basis functions and gradients in a cache
Definition: localbasiscache.hh:18
LocalBasisCache()
constructor
Definition: localbasiscache.hh:43
const std::vector< RangeType > & evaluateFunction(const DomainType &position, const LocalBasisType &localbasis) const
evaluate basis functions at a point
Definition: localbasiscache.hh:47
const std::vector< JacobianType > & evaluateJacobian(const DomainType &position, const LocalBasisType &localbasis) const
evaluate Jacobians at a point
Definition: localbasiscache.hh:59
A few common exception classes.
Dune namespace.
Definition: alignedallocator.hh:14
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)