Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d6/dbe/Vector_8h_source.html
DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
Vector.h
1// Copyright (C) 2020 Garth N. Wells
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#include "utils.h"
10#include <cmath>
11#include <complex>
12#include <dolfinx/common/IndexMap.h>
13#include <dolfinx/common/Scatterer.h>
14#include <dolfinx/common/types.h>
15#include <limits>
16#include <memory>
17#include <numeric>
18#include <span>
19#include <type_traits>
20#include <vector>
21
22namespace dolfinx::la
23{
24
29template <typename T, typename Container = std::vector<T>>
30class Vector
31{
32 static_assert(std::is_same_v<typename Container::value_type, T>);
33
34public:
36 using value_type = T;
37
39 using container_type = Container;
40
41 static_assert(std::is_same_v<value_type, typename container_type::value_type>,
42 "Scalar type and container value type must be the same.");
43
47 Vector(std::shared_ptr<const common::IndexMap> map, int bs)
48 : _map(map), _scatterer(std::make_shared<common::Scatterer<>>(*_map, bs)),
49 _bs(bs), _buffer_local(_scatterer->local_buffer_size()),
50 _buffer_remote(_scatterer->remote_buffer_size()),
51 _x(bs * (map->size_local() + map->num_ghosts()))
52 {
53 }
54
56 Vector(const Vector& x)
57 : _map(x._map), _scatterer(x._scatterer), _bs(x._bs),
58 _request(1, MPI_REQUEST_NULL), _buffer_local(x._buffer_local),
59 _buffer_remote(x._buffer_remote), _x(x._x)
60 {
61 }
62
65 : _map(std::move(x._map)), _scatterer(std::move(x._scatterer)),
66 _bs(std::move(x._bs)),
67 _request(std::exchange(x._request, {MPI_REQUEST_NULL})),
68 _buffer_local(std::move(x._buffer_local)),
69 _buffer_remote(std::move(x._buffer_remote)), _x(std::move(x._x))
70 {
71 }
72
73 // Assignment operator (disabled)
74 Vector& operator=(const Vector& x) = delete;
75
77 Vector& operator=(Vector&& x) = default;
78
81 void set(value_type v) { std::fill(_x.begin(), _x.end(), v); }
82
86 {
87 const std::int32_t local_size = _bs * _map->size_local();
88 std::span<const value_type> x_local(_x.data(), local_size);
89
90 auto pack = [](auto&& in, auto&& idx, auto&& out)
91 {
92 for (std::size_t i = 0; i < idx.size(); ++i)
93 out[i] = in[idx[i]];
94 };
95 pack(x_local, _scatterer->local_indices(), _buffer_local);
96
97 _scatterer->scatter_fwd_begin(std::span<const value_type>(_buffer_local),
98 std::span<value_type>(_buffer_remote),
99 std::span<MPI_Request>(_request));
100 }
101
105 {
106 const std::int32_t local_size = _bs * _map->size_local();
107 const std::int32_t num_ghosts = _bs * _map->num_ghosts();
108 std::span<value_type> x_remote(_x.data() + local_size, num_ghosts);
109 _scatterer->scatter_fwd_end(std::span<MPI_Request>(_request));
110
111 auto unpack = [](auto&& in, auto&& idx, auto&& out, auto op)
112 {
113 for (std::size_t i = 0; i < idx.size(); ++i)
114 out[idx[i]] = op(out[idx[i]], in[i]);
115 };
116
117 unpack(_buffer_remote, _scatterer->remote_indices(), x_remote,
118 [](auto /*a*/, auto b) { return b; });
119 }
120
124 {
125 this->scatter_fwd_begin();
126 this->scatter_fwd_end();
127 }
128
132 {
133 const std::int32_t local_size = _bs * _map->size_local();
134 const std::int32_t num_ghosts = _bs * _map->num_ghosts();
135 std::span<value_type> x_remote(_x.data() + local_size, num_ghosts);
136
137 auto pack = [](auto&& in, auto&& idx, auto&& out)
138 {
139 for (std::size_t i = 0; i < idx.size(); ++i)
140 out[i] = in[idx[i]];
141 };
142 pack(x_remote, _scatterer->remote_indices(), _buffer_remote);
143
144 _scatterer->scatter_rev_begin(std::span<const value_type>(_buffer_remote),
145 std::span<value_type>(_buffer_local),
146 _request);
147 }
148
155 template <class BinaryOperation>
156 void scatter_rev_end(BinaryOperation op)
157 {
158 const std::int32_t local_size = _bs * _map->size_local();
159 std::span<value_type> x_local(_x.data(), local_size);
160 _scatterer->scatter_rev_end(_request);
161
162 auto unpack = [](auto&& in, auto&& idx, auto&& out, auto op)
163 {
164 for (std::size_t i = 0; i < idx.size(); ++i)
165 out[idx[i]] = op(out[idx[i]], in[i]);
166 };
167 unpack(_buffer_local, _scatterer->local_indices(), x_local, op);
168 }
169
175 template <class BinaryOperation>
176 void scatter_rev(BinaryOperation op)
177 {
178 this->scatter_rev_begin();
179 this->scatter_rev_end(op);
180 }
181
183 std::shared_ptr<const common::IndexMap> index_map() const { return _map; }
184
186 constexpr int bs() const { return _bs; }
187
189 std::span<const value_type> array() const
190 {
191 return std::span<const value_type>(_x);
192 }
193
195 std::span<value_type> mutable_array() { return std::span(_x); }
196
197private:
198 // Map describing the data layout
199 std::shared_ptr<const common::IndexMap> _map;
200
201 // Scatter for managing MPI communication
202 std::shared_ptr<const common::Scatterer<>> _scatterer;
203
204 // Block size
205 int _bs;
206
207 // MPI request handle
208 std::vector<MPI_Request> _request = {MPI_REQUEST_NULL};
209
210 // Buffers for ghost scatters
211 container_type _buffer_local, _buffer_remote;
212
213 // Vector data
215};
216
223template <class V>
224auto inner_product(const V& a, const V& b)
225{
226 using T = typename V::value_type;
227 const std::int32_t local_size = a.bs() * a.index_map()->size_local();
228 if (local_size != b.bs() * b.index_map()->size_local())
229 throw std::runtime_error("Incompatible vector sizes");
230 std::span<const T> x_a = a.array().subspan(0, local_size);
231 std::span<const T> x_b = b.array().subspan(0, local_size);
232
233 const T local = std::transform_reduce(
234 x_a.begin(), x_a.end(), x_b.begin(), static_cast<T>(0), std::plus{},
235 [](T a, T b) -> T
236 {
237 if constexpr (std::is_same<T, std::complex<double>>::value
238 or std::is_same<T, std::complex<float>>::value)
239 {
240 return std::conj(a) * b;
241 }
242 else
243 return a * b;
244 });
245
246 T result;
247 MPI_Allreduce(&local, &result, 1, dolfinx::MPI::mpi_type<T>(), MPI_SUM,
248 a.index_map()->comm());
249 return result;
250}
251
254template <class V>
255auto squared_norm(const V& a)
256{
257 using T = typename V::value_type;
258 T result = inner_product(a, a);
259 return std::real(result);
260}
261
266template <class V>
267auto norm(const V& x, Norm type = Norm::l2)
268{
269 using T = typename V::value_type;
270 switch (type)
271 {
272 case Norm::l1:
273 {
274 std::int32_t size_local = x.bs() * x.index_map()->size_local();
275 std::span<const T> data = x.array().subspan(0, size_local);
276 using U = typename dolfinx::scalar_value_type_t<T>;
277 U local_l1
278 = std::accumulate(data.begin(), data.end(), U(0),
279 [](auto norm, auto x) { return norm + std::abs(x); });
280 U l1(0);
281 MPI_Allreduce(&local_l1, &l1, 1, MPI::mpi_type<U>(), MPI_SUM,
282 x.index_map()->comm());
283 return l1;
284 }
285 case Norm::l2:
286 return std::sqrt(squared_norm(x));
287 case Norm::linf:
288 {
289 std::int32_t size_local = x.bs() * x.index_map()->size_local();
290 std::span<const T> data = x.array().subspan(0, size_local);
291 auto max_pos = std::max_element(data.begin(), data.end(), [](T a, T b)
292 { return std::norm(a) < std::norm(b); });
293 auto local_linf = std::abs(*max_pos);
294 decltype(local_linf) linf = 0;
295 MPI_Allreduce(&local_linf, &linf, 1, MPI::mpi_type<decltype(linf)>(),
296 MPI_MAX, x.index_map()->comm());
297 return linf;
298 }
299 default:
300 throw std::runtime_error("Norm type not supported");
301 }
302}
303
309template <class V>
310void orthonormalize(std::vector<std::reference_wrapper<V>> basis)
311{
312 using T = typename V::value_type;
313 using U = typename dolfinx::scalar_value_type_t<T>;
314
315 // Loop over each vector in basis
316 for (std::size_t i = 0; i < basis.size(); ++i)
317 {
318 // Orthogonalize vector i with respect to previously orthonormalized
319 // vectors
320 V& bi = basis[i].get();
321 for (std::size_t j = 0; j < i; ++j)
322 {
323 const V& bj = basis[j].get();
324
325 // basis_i <- basis_i - dot_ij basis_j
326 auto dot_ij = inner_product(bi, bj);
327 std::transform(bj.array().begin(), bj.array().end(), bi.array().begin(),
328 bi.mutable_array().begin(),
329 [dot_ij](auto xj, auto xi) { return xi - dot_ij * xj; });
330 }
331
332 // Normalise basis function
333 auto norm = la::norm(bi, la::Norm::l2);
334 if (norm * norm < std::numeric_limits<U>::epsilon())
335 {
336 throw std::runtime_error(
337 "Linear dependency detected. Cannot orthogonalize.");
338 }
339 std::transform(bi.array().begin(), bi.array().end(),
340 bi.mutable_array().begin(),
341 [norm](auto x) { return x / norm; });
342 }
343}
344
348template <class V>
349bool is_orthonormal(std::vector<std::reference_wrapper<const V>> basis)
350{
351 using T = typename V::value_type;
352 using U = typename dolfinx::scalar_value_type_t<T>;
353
354 // auto tol = std::sqrt(T(std::numeric_limits<U>::epsilon()));
355 for (std::size_t i = 0; i < basis.size(); i++)
356 {
357 for (std::size_t j = i; j < basis.size(); j++)
358 {
359 T delta_ij = (i == j) ? T(1) : T(0);
360 auto dot_ij = inner_product(basis[i].get(), basis[j].get());
361 if (std::norm(delta_ij - dot_ij) > std::numeric_limits<U>::epsilon())
362 return false;
363 }
364 }
365
366 return true;
367}
368
369} // namespace dolfinx::la
Definition Vector.h:31
void scatter_fwd_begin()
Definition Vector.h:85
constexpr int bs() const
Get block size.
Definition Vector.h:186
Vector(std::shared_ptr< const common::IndexMap > map, int bs)
Definition Vector.h:47
Vector & operator=(Vector &&x)=default
Move Assignment operator.
void scatter_fwd()
Definition Vector.h:123
std::span< const value_type > array() const
Get local part of the vector (const version)
Definition Vector.h:189
std::shared_ptr< const common::IndexMap > index_map() const
Get IndexMap.
Definition Vector.h:183
Container container_type
Container type.
Definition Vector.h:39
Vector(Vector &&x)
Move constructor.
Definition Vector.h:64
Vector(const Vector &x)
Copy constructor.
Definition Vector.h:56
void scatter_fwd_end()
Definition Vector.h:104
void scatter_rev_begin()
Definition Vector.h:131
void set(value_type v)
Definition Vector.h:81
std::span< value_type > mutable_array()
Get local part of the vector.
Definition Vector.h:195
void scatter_rev_end(BinaryOperation op)
Definition Vector.h:156
void scatter_rev(BinaryOperation op)
Definition Vector.h:176
T value_type
Scalar type.
Definition Vector.h:36
Linear algebra interface.
Definition sparsitybuild.h:15
void orthonormalize(std::vector< std::reference_wrapper< V > > basis)
Definition Vector.h:310
Norm
Norm types.
Definition utils.h:17
auto squared_norm(const V &a)
Definition Vector.h:255
bool is_orthonormal(std::vector< std::reference_wrapper< const V > > basis)
Test if basis is orthonormal.
Definition Vector.h:349
auto norm(const V &x, Norm type=Norm::l2)
Definition Vector.h:267
auto inner_product(const V &a, const V &b)
Definition Vector.h:224