11 #include <dolfinx/common/IndexMap.h>
12 #include <dolfinx/common/Scatterer.h>
24 template <
typename T,
class Allocator = std::allocator<T>>
35 Vector(
const std::shared_ptr<const common::IndexMap>&
map,
int bs,
36 const Allocator& alloc = Allocator())
37 : _map(
map), _scatterer(std::make_shared<common::Scatterer>(*_map,
bs)),
38 _bs(
bs), _buffer_local(_scatterer->local_buffer_size(), alloc),
39 _buffer_remote(_scatterer->remote_buffer_size(), alloc),
40 _x(
bs * (
map->size_local() +
map->num_ghosts()), alloc)
46 : _map(x._map), _scatterer(x._scatterer), _bs(x._bs),
47 _request(MPI_REQUEST_NULL), _buffer_local(x._buffer_local),
48 _buffer_remote(x._buffer_remote), _x(x._x)
54 : _map(std::move(x._map)), _scatterer(std::move(x._scatterer)),
55 _bs(std::move(x._bs)),
56 _request(std::exchange(x._request, MPI_REQUEST_NULL)),
57 _buffer_local(std::move(x._buffer_local)),
58 _buffer_remote(std::move(x._buffer_remote)), _x(std::move(x._x))
70 void set(T v) { std::fill(_x.begin(), _x.end(), v); }
76 const std::int32_t local_size = _bs * _map->size_local();
77 std::span<const T> x_local(_x.data(), local_size);
79 auto pack = [](
const auto& in,
const auto& idx,
auto& out)
81 for (std::size_t i = 0; i < idx.size(); ++i)
84 pack(x_local, _scatterer->local_indices(), _buffer_local);
86 _scatterer->scatter_fwd_begin(std::span<const T>(_buffer_local),
87 std::span<T>(_buffer_remote), _request);
94 const std::int32_t local_size = _bs * _map->size_local();
95 const std::int32_t num_ghosts = _bs * _map->num_ghosts();
96 std::span<T> x_remote(_x.data() + local_size, num_ghosts);
97 _scatterer->scatter_fwd_end(_request);
99 auto unpack = [](
const auto& in,
const auto& idx,
auto& out,
auto op)
101 for (std::size_t i = 0; i < idx.size(); ++i)
102 out[idx[i]] = op(out[idx[i]], in[i]);
105 unpack(_buffer_remote, _scatterer->remote_indices(), x_remote,
106 [](
auto ,
auto b) { return b; });
121 const std::int32_t local_size = _bs * _map->size_local();
122 const std::int32_t num_ghosts = _bs * _map->num_ghosts();
123 std::span<T> x_remote(_x.data() + local_size, num_ghosts);
125 auto pack = [](
const auto& in,
const auto& idx,
auto& out)
127 for (std::size_t i = 0; i < idx.size(); ++i)
130 pack(x_remote, _scatterer->remote_indices(), _buffer_remote);
132 _scatterer->scatter_rev_begin(std::span<const T>(_buffer_remote),
133 std::span<T>(_buffer_local), _request);
142 template <
class BinaryOperation>
145 const std::int32_t local_size = _bs * _map->size_local();
146 std::span<T> x_local(_x.data(), local_size);
147 _scatterer->scatter_rev_end(_request);
149 auto unpack = [](
const auto& in,
const auto& idx,
auto& out,
auto op)
151 for (std::size_t i = 0; i < idx.size(); ++i)
152 out[idx[i]] = op(out[idx[i]], in[i]);
154 unpack(_buffer_local, _scatterer->local_indices(), x_local, op);
162 template <
class BinaryOperation>
170 std::shared_ptr<const common::IndexMap>
map()
const {
return _map; }
173 constexpr
int bs()
const {
return _bs; }
176 std::span<const T>
array()
const {
return std::span<const T>(_x); }
186 std::shared_ptr<const common::IndexMap> _map;
189 std::shared_ptr<const common::Scatterer> _scatterer;
195 MPI_Request _request = MPI_REQUEST_NULL;
198 std::vector<T, Allocator> _buffer_local, _buffer_remote;
201 std::vector<T, Allocator> _x;
210 template <
typename T,
class Allocator>
213 const std::int32_t local_size = a.
bs() * a.
map()->size_local();
214 if (local_size != b.
bs() * b.
map()->size_local())
215 throw std::runtime_error(
"Incompatible vector sizes");
216 std::span<const T> x_a = a.
array().subspan(0, local_size);
217 std::span<const T> x_b = b.
array().subspan(0, local_size);
219 const T local = std::transform_reduce(
220 x_a.begin(), x_a.end(), x_b.begin(),
static_cast<T
>(0), std::plus{},
223 if constexpr (std::is_same<T, std::complex<double>>::value
224 or std::is_same<T, std::complex<float>>::value)
226 return std::conj(a) * b;
233 MPI_Allreduce(&local, &result, 1, dolfinx::MPI::mpi_type<T>(), MPI_SUM,
240 template <
typename T,
class Allocator>
244 return std::real(result);
251 template <
typename T,
class Allocator>
260 const std::int32_t size_local = a.
bs() * a.
map()->size_local();
261 std::span<const T> x_a = a.
array().subspan(0, size_local);
262 auto max_pos = std::max_element(x_a.begin(), x_a.end(),
264 { return std::norm(a) < std::norm(b); });
265 auto local_linf = std::abs(*max_pos);
266 decltype(local_linf) linf = 0;
267 MPI_Allreduce(&local_linf, &linf, 1, MPI::mpi_type<decltype(linf)>(),
268 MPI_MAX, a.map()->comm());
272 throw std::runtime_error(
"Norm type not supported");
281 template <
typename T,
typename U>
285 for (std::size_t i = 0; i < basis.size(); ++i)
289 for (std::size_t j = 0; j < i; ++j)
293 std::transform(basis[j].array().begin(), basis[j].array().end(),
294 basis[i].array().begin(), basis[i].mutable_array().begin(),
295 [dot_ij](
auto xj,
auto xi) {
return xi - dot_ij * xj; });
299 double norm = la::norm(basis[i], la::Norm::l2);
302 throw std::runtime_error(
303 "Linear dependency detected. Cannot orthogonalize.");
305 std::transform(basis[i].array().begin(), basis[i].array().end(),
306 basis[i].mutable_array().begin(),
307 [
norm](
auto x) {
return x /
norm; });
315 template <
typename T,
typename U>
317 double tol = 1.0e-10)
319 for (std::size_t i = 0; i < basis.size(); i++)
321 for (std::size_t j = i; j < basis.size(); j++)
323 const double delta_ij = (i == j) ? 1.0 : 0.0;
325 if (std::abs(delta_ij - dot_ij) > tol)
Distributed vector.
Definition: Vector.h:26
std::shared_ptr< const common::IndexMap > map() const
Get IndexMap.
Definition: Vector.h:170
void scatter_fwd_begin()
Begin scatter of local data from owner to ghosts on other ranks.
Definition: Vector.h:74
Allocator allocator_type
The allocator type.
Definition: Vector.h:32
std::span< const T > array() const
Get local part of the vector (const version)
Definition: Vector.h:176
Vector(const std::shared_ptr< const common::IndexMap > &map, int bs, const Allocator &alloc=Allocator())
Create a distributed vector.
Definition: Vector.h:35
constexpr int bs() const
Get block size.
Definition: Vector.h:173
void set(T v)
Set all entries (including ghosts)
Definition: Vector.h:70
Vector & operator=(Vector &&x)=default
Move Assignment operator.
void scatter_fwd()
Scatter local data to ghost positions on other ranks.
Definition: Vector.h:111
Vector(Vector &&x)
Move constructor.
Definition: Vector.h:53
Vector(const Vector &x)
Copy constructor.
Definition: Vector.h:45
void scatter_fwd_end()
End scatter of local data from owner to ghosts on other ranks.
Definition: Vector.h:92
void scatter_rev_begin()
Start scatter of ghost data to owner.
Definition: Vector.h:119
constexpr allocator_type allocator() const
Get the allocator associated with the container.
Definition: Vector.h:182
std::span< T > mutable_array()
Get local part of the vector.
Definition: Vector.h:179
void scatter_rev_end(BinaryOperation op)
End scatter of ghost data to owner. This process may receive data from more than one process,...
Definition: Vector.h:143
void scatter_rev(BinaryOperation op)
Scatter ghost data to owner. This process may receive data from more than one process,...
Definition: Vector.h:163
T value_type
The value type.
Definition: Vector.h:29
Linear algebra interface.
Definition: sparsitybuild.h:15
Norm
Norm types.
Definition: utils.h:13
auto norm(const Vector< T, Allocator > &a, Norm type=Norm::l2)
Compute the norm of the vector.
Definition: Vector.h:252
T inner_product(const Vector< T, Allocator > &a, const Vector< T, Allocator > &b)
Compute the inner product of two vectors. The two vectors must have the same parallel layout.
Definition: Vector.h:211
void orthonormalize(const std::span< Vector< T, U >> &basis, double tol=1.0e-10)
Orthonormalize a set of vectors.
Definition: Vector.h:282
auto squared_norm(const Vector< T, Allocator > &a)
Compute the squared L2 norm of vector.
Definition: Vector.h:241
bool is_orthonormal(const std::span< const Vector< T, U >> &basis, double tol=1.0e-10)
Test if basis is orthonormal.
Definition: Vector.h:316