DUNE-FEM (unstable)

pqkfactory.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5#ifndef DUNE_LOCALFUNCTIONS_LAGRANGE_PQKFACTORY_HH
6#define DUNE_LOCALFUNCTIONS_LAGRANGE_PQKFACTORY_HH
7
8#include <map>
9
10#include <dune/geometry/type.hh>
11
12#include <dune/localfunctions/common/virtualinterface.hh>
13#include <dune/localfunctions/common/virtualwrappers.hh>
14
15#include <dune/localfunctions/lagrange/lagrangecube.hh>
16#include <dune/localfunctions/lagrange/lagrangeprism.hh>
17#include <dune/localfunctions/lagrange/lagrangepyramid.hh>
18#include <dune/localfunctions/lagrange/lagrangesimplex.hh>
19#include <dune/localfunctions/lagrange/p0.hh>
20
21namespace Dune
22{
23
28 template<class D, class R, int d, int k>
30 {
32
35 {
36 return nullptr;
37 }
38 };
39
44 template<class D, class R, int k>
46 {
52
55 {
56 if ((gt.isPrism())and (k==1))
58 if ((gt.isPrism())and (k==2))
60 if ((gt.isPyramid())and (k==1))
62 if ((gt.isPyramid())and (k==2))
64 return nullptr;
65 }
66 };
67
68
72 template<class D, class R, int dim, int k>
74 {
80
81
84 {
85 if (k==0)
87
88 if (gt.isSimplex())
90
91 if (gt.isCube())
93
95 }
96 };
97
98
99
110 template<class D, class R, int dim, int k>
111 class
112 [[deprecated("Use LagrangeLocalFiniteElementCache<D,R,dim,k> from lagrangelfecache.hh. This will be removed after release 2.10.")]]
114 {
115 protected:
118 typedef typename std::map<GeometryType,FE*> FEMap;
119
120 public:
123
126
129 {
130 typename FEMap::iterator it = other.cache_.begin();
131 typename FEMap::iterator end = other.cache_.end();
132 for(; it!=end; ++it)
133 cache_[it->first] = (it->second)->clone();
134 }
135
137 {
138 typename FEMap::iterator it = cache_.begin();
139 typename FEMap::iterator end = cache_.end();
140 for(; it!=end; ++it)
141 delete it->second;
142 }
143
146 {
147 typename FEMap::const_iterator it = cache_.find(gt);
148 if (it==cache_.end())
149 {
151 if (fe==0)
152 DUNE_THROW(Dune::NotImplemented,"No Pk/Qk like local finite element available for geometry type " << gt << " and order " << k);
153
154 cache_[gt] = fe;
155 return *fe;
156 }
157 return *(it->second);
158 }
159
160 protected:
161 mutable FEMap cache_;
162
163 };
164
165}
166
167#endif
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
Lagrange finite element for cubes with arbitrary compile-time dimension and polynomial order.
Definition: lagrangecube.hh:709
Lagrange finite element for 3d prisms with arbitrary compile-time polynomial order.
Definition: lagrangeprism.hh:649
Lagrange finite element for 3d pyramids with compile-time polynomial order.
Definition: lagrangepyramid.hh:809
Lagrange finite element for simplices with arbitrary compile-time dimension and polynomial order.
Definition: lagrangesimplex.hh:836
class for wrapping a finite element using the virtual interface
Definition: virtualwrappers.hh:240
virtual base class for local finite elements with functions
Definition: virtualinterface.hh:225
Default exception for dummy implementations.
Definition: exceptions.hh:263
The local p0 finite element on all types of reference elements.
Definition: p0.hh:25
A cache that stores all available Pk/Qk like local finite elements for the given dimension and order.
Definition: pqkfactory.hh:114
PQkLocalFiniteElementCache()
Default constructor.
Definition: pqkfactory.hh:125
FE FiniteElementType
Type of the finite elements stored in this cache.
Definition: pqkfactory.hh:122
const FiniteElementType & get(const GeometryType &gt) const
Get local finite element for given GeometryType.
Definition: pqkfactory.hh:145
PQkLocalFiniteElementCache(const PQkLocalFiniteElementCache &other)
Copy constructor.
Definition: pqkfactory.hh:128
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:158
Dune namespace.
Definition: alignedallocator.hh:13
static LocalFiniteElementVirtualInterface< T > * create(const GeometryType &gt)
create finite element for given GeometryType
Definition: pqkfactory.hh:54
Factory that only creates dimension specific local finite elements.
Definition: pqkfactory.hh:30
static LocalFiniteElementVirtualInterface< T > * create(const GeometryType &)
create finite element for given GeometryType
Definition: pqkfactory.hh:34
Factory to create any kind of Pk/Qk like element wrapped for the virtual interface.
Definition: pqkfactory.hh:74
static FiniteElementType * create(const GeometryType &gt)
create finite element for given GeometryType
Definition: pqkfactory.hh:83
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)