3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
reservedblockvector.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_RESERVED_BLOCK_VECTOR_HH
25#define DUMUX_RESERVED_BLOCK_VECTOR_HH
26
27#include <algorithm>
28#include <dune/common/reservedvector.hh>
29
30namespace Dumux {
31
36template<class BlockType, int capacity>
37class ReservedBlockVector : public Dune::ReservedVector<BlockType, capacity>
38{
39 using Base = Dune::ReservedVector<BlockType, capacity>;
40public:
41
42 using size_type = typename Base::size_type;
43 using value_type = BlockType;
44
45 using Base::Base;
46
47 explicit ReservedBlockVector() : Base() {}
48 explicit ReservedBlockVector(size_type size) : Base() { this->resize(size); }
49
52
55
57
59 ReservedBlockVector& operator= (const typename BlockType::field_type& v)
60 {
61 std::fill(this->begin(), this->end(), v);
62 return *this;
63 }
64
67 {
68 for (size_type i = 0; i < this->size(); ++i)
69 (*this)[i] += other[i];
70 return *this;
71 }
72
75 {
76 for (size_type i = 0; i < this->size(); ++i)
77 (*this)[i] -= other[i];
78 return *this;
79 }
80
82 ReservedBlockVector& operator/= (const typename BlockType::field_type& v)
83 {
84 for (size_type i = 0; i < this->size(); ++i)
85 (*this)[i] /= v;
86 return *this;
87 }
88
90 ReservedBlockVector& operator*= (const typename BlockType::field_type& v)
91 {
92 for (size_type i = 0; i < this->size(); ++i)
93 (*this)[i] *= v;
94 return *this;
95 }
96};
97
98} // end namespace Dumux
99
100#endif
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
void resize(const Vector &v, std::size_t size)
Definition: test_isvalid.cc:26
A arithmetic block vector type based on DUNE's reserved vector.
Definition: reservedblockvector.hh:38
ReservedBlockVector(const ReservedBlockVector &)=default
ReservedBlockVector & operator=(const ReservedBlockVector &)=default
ReservedBlockVector & operator-=(const ReservedBlockVector &other)
vector space subtraction
Definition: reservedblockvector.hh:74
ReservedBlockVector & operator+=(const ReservedBlockVector &other)
vector space addition
Definition: reservedblockvector.hh:66
ReservedBlockVector & operator/=(const typename BlockType::field_type &v)
division by scalar
Definition: reservedblockvector.hh:82
typename Base::size_type size_type
Definition: reservedblockvector.hh:42
ReservedBlockVector(ReservedBlockVector &&)=default
ReservedBlockVector()
Definition: reservedblockvector.hh:47
BlockType value_type
Definition: reservedblockvector.hh:43
ReservedBlockVector(size_type size)
Definition: reservedblockvector.hh:48
ReservedBlockVector & operator*=(const typename BlockType::field_type &v)
multiplication by scalar
Definition: reservedblockvector.hh:90