Dune Core Modules (2.10.0)

interpolationhelper.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 GENERIC_INTERPOLATIONHELPER_HH
6#define GENERIC_INTERPOLATIONHELPER_HH
7
8#include <vector>
9
12#include <dune/localfunctions/utility/field.hh>
13
14namespace Dune
15{
16 // A small helper class to avoid having to
17 // write the interpolation twice (once for function
18 // and once for a basis)
19 template< class F, unsigned int dimension >
20 struct InterpolationHelper
21 {
22 template <class Func,class Container, bool type>
23 struct Helper;
24 };
25 template <class F,unsigned int d>
26 template <class Func,class Vector>
27 struct InterpolationHelper<F,d>::Helper<Func,Vector,true>
28 // Func is of Function type
29 {
30 typedef std::vector< Dune::FieldVector<F,d> > Result;
31 Helper(const Func & func, Vector &vec)
32 : func_(func),
33 vec_(vec),
34 tmp_(1)
35 {}
36 const typename Vector::value_type &operator()(unsigned int row,unsigned int col)
37 {
38 return vec_[row];
39 }
40 template <class Fy>
41 void set(unsigned int row,unsigned int col,
42 const Fy &val)
43 {
44 assert(col==0);
45 assert(row<vec_.size());
46 field_cast( val, vec_[row] );
47 }
48 template <class Fy>
49 void add(unsigned int row,unsigned int col,
50 const Fy &val)
51 {
52 assert(col==0);
53 assert(row<vec_.size());
54 vec_[row] += field_cast<typename Vector::value_type>(val);
55 }
56 template <class DomainVector>
57 const Result &evaluate(const DomainVector &x) const
58 {
59 field_cast(func_(x), tmp_[0] );
60 return tmp_;
61 }
62 unsigned int size() const
63 {
64 return 1;
65 }
66 const Func &func_;
67 Vector &vec_;
68 mutable Result tmp_;
69 };
70 template <class F,unsigned int d>
71 template <class Basis,class Matrix>
72 struct InterpolationHelper<F,d>::Helper<Basis,Matrix,false>
73 // Func is of Basis type
74 {
75 typedef std::vector< Dune::FieldVector<F,d> > Result;
76 Helper(const Basis & basis, Matrix &matrix)
77 : basis_(basis),
78 matrix_(matrix),
79 tmp_(basis.size()) {}
80 const F &operator()(unsigned int row,unsigned int col) const
81 {
82 return matrix_[row][col];
83 }
84 F &operator()(unsigned int row,unsigned int col)
85 {
86 return matrix_[row][col];
87 }
88 template <class Fy>
89 void set(unsigned int row,unsigned int col,
90 const Fy &val)
91 {
92 assert(col<matrix_.cols());
93 assert(row<matrix_.rows());
94 field_cast(val,matrix_[row][col]);
95 }
96 template <class Fy>
97 void add(unsigned int row,unsigned int col,
98 const Fy &val)
99 {
100 assert(col<matrix_.cols());
101 assert(row<matrix_.rows());
102 matrix_[row][col] += val;
103 }
104 template <class DomainVector>
105 const Result &evaluate(const DomainVector &x) const
106 {
107 basis_.template evaluate<0>(x,tmp_);
108 return tmp_;
109 }
110 unsigned int size() const
111 {
112 return basis_.size();
113 }
114
115 const Basis &basis_;
116 Matrix &matrix_;
117 mutable Result tmp_;
118 };
119}
120#endif // GENERIC_INTERPOLATIONHELPER_HH
Infrastructure for concepts.
Implements a vector constructed from a given type representing a field and a compile-time given size.
Dune namespace.
Definition: alignedallocator.hh:13
void field_cast(const F1 &f1, F2 &f2)
a helper class to cast from one field to another
Definition: field.hh:159
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 28, 23:30, 2024)