13#include <dolfinx/common/IndexMap.h>
14#include <dolfinx/common/Scatterer.h>
15#include <dolfinx/common/types.h>
26template <
class F,
class Container,
class ScatterContainer>
28 f(idx.cbegin(), idx.cend(), x.cbegin(), x.begin());
32template <
class GetPtr,
class U,
class T>
34 { f(x) } -> std::same_as<T*>;
35} and
requires(GetPtr f,
const U x) {
36 { f(x) } -> std::same_as<const T*>;
47template <
typename T,
typename Container = std::vector<T>,
48 typename ScatterContainer = std::vector<std::
int32_t>>
51 static_assert(std::is_same_v<typename Container::value_type, T>);
60 return [](
typename ScatterContainer::const_iterator idx_first,
61 typename ScatterContainer::const_iterator idx_last,
62 const auto in_first,
auto out_first)
65 std::transform(idx_first, idx_last, out_first,
66 [in_first](
auto p) {
return *std::next(in_first, p); });
76 return [](
typename ScatterContainer::const_iterator idx_first,
77 typename ScatterContainer::const_iterator idx_last,
78 const auto in_first,
auto out_first)
81 for (
typename ScatterContainer::const_iterator idx = idx_first;
82 idx != idx_last; ++idx)
84 std::size_t d = std::distance(idx_first, idx);
85 *std::next(out_first, *idx) = *std::next(in_first, d);
95 template <
typename BinaryOp>
96 auto get_unpack_op(BinaryOp op)
98 return [op](
typename ScatterContainer::const_iterator idx_first,
99 typename ScatterContainer::const_iterator idx_last,
100 const auto in_first,
auto out_first)
103 for (
typename ScatterContainer::const_iterator idx = idx_first;
104 idx != idx_last; ++idx)
106 std::size_t d = std::distance(idx_first, idx);
107 auto& out = *std::next(out_first, *idx);
108 out = op(out, *std::next(in_first, d));
120 static_assert(std::is_same_v<value_type, typename container_type::value_type>,
121 "Scalar type and container value type must be the same.");
128 Vector(std::shared_ptr<const common::IndexMap> map,
int bs)
129 : _map(map), _bs(
bs), _x(
bs * (map->size_local() + map->num_ghosts())),
131 std::make_shared<
common::Scatterer<ScatterContainer>>(*_map,
bs)),
132 _buffer_local(_scatterer->local_indices().size()),
133 _buffer_remote(_scatterer->remote_indices().size())
153 [[deprecated(
"Use std::ranges::fill(u.array(), v) instead.")]]
void
156 std::ranges::fill(_x, v);
173 template <
typename U,
typename GetPtr>
178 pack(_scatterer->local_indices().begin(), _scatterer->local_indices().end(),
179 _x.begin(), _buffer_local.begin());
180 _scatterer->scatter_fwd_begin(get_ptr(_buffer_local),
181 get_ptr(_buffer_remote), _request);
193 requires requires(Container c) {
194 { c.data() } -> std::same_as<T*>;
211 template <
typename U>
212 requires VectorPackKernel<U, container_type, ScatterContainer>
215 _scatterer->scatter_end(_request);
216 unpack(_scatterer->remote_indices().begin(),
217 _scatterer->remote_indices().end(), _buffer_remote.begin(),
218 std::next(_x.begin(), _bs * _map->size_local()));
230 requires requires(Container c) {
231 { c.data() } -> std::same_as<T*>;
247 requires requires(Container c) {
248 { c.data() } -> std::same_as<T*>;
268 template <
typename U,
typename GetPtr>
269 requires VectorPackKernel<U, container_type, ScatterContainer>
270 and GetPtrConcept<GetPtr, Container, T>
273 std::int32_t local_size = _bs * _map->size_local();
274 pack(_scatterer->remote_indices().begin(),
275 _scatterer->remote_indices().end(), std::next(_x.begin(), local_size),
276 _buffer_remote.begin());
277 _scatterer->scatter_rev_begin(get_ptr(_buffer_remote),
278 get_ptr(_buffer_local), _request);
290 requires requires(Container c) {
291 { c.data() } -> std::same_as<T*>;
304 template <
typename U>
305 requires VectorPackKernel<U, container_type, ScatterContainer>
308 _scatterer->scatter_end(_request);
309 unpack(_scatterer->local_indices().begin(),
310 _scatterer->local_indices().end(), _buffer_local.begin(),
325 template <
class BinaryOperation>
326 requires requires(Container c) {
327 { c.data() } -> std::same_as<T*>;
336 std::shared_ptr<const common::IndexMap>
index_map()
const {
return _map; }
339 constexpr int bs()
const {
return _bs; }
360 std::shared_ptr<const common::IndexMap> _map;
369 std::shared_ptr<const common::Scatterer<ScatterContainer>> _scatterer;
372 MPI_Request _request = MPI_REQUEST_NULL;
391 using T =
typename V::value_type;
392 const std::int32_t local_size = a.bs() * a.index_map()->size_local();
393 if (local_size != b.bs() * b.index_map()->size_local())
394 throw std::runtime_error(
"Incompatible vector sizes");
396 const T local = std::transform_reduce(
397 a.array().begin(), std::next(a.array().begin(), local_size),
398 b.array().begin(),
static_cast<T
>(0), std::plus{},
401 if constexpr (std::is_same<T, std::complex<double>>::value
402 or std::is_same<T, std::complex<float>>::value)
404 return std::conj(a) * b;
412 a.index_map()->comm());
422 using T =
typename V::value_type;
424 return std::real(result);
436 using T =
typename V::value_type;
441 std::int32_t size_local = x.bs() * x.index_map()->size_local();
442 using U =
typename dolfinx::scalar_value_t<T>;
443 U local_l1 = std::accumulate(
444 x.array().begin(), std::next(x.array().begin(), size_local), U(0),
445 [](
auto norm,
auto x) {
return norm + std::abs(x); });
448 x.index_map()->comm());
455 std::int32_t size_local = x.bs() * x.index_map()->size_local();
456 auto max_pos = std::max_element(
457 x.array().begin(), std::next(x.array().begin(), size_local),
458 [](T a, T b) { return std::norm(a) < std::norm(b); });
459 auto local_linf = std::abs(*max_pos);
460 decltype(local_linf) linf = 0;
461 MPI_Allreduce(&local_linf, &linf, 1,
MPI::mpi_t<
decltype(linf)>, MPI_MAX,
462 x.index_map()->comm());
466 throw std::runtime_error(
"Norm type not supported");
479 using T =
typename V::value_type;
480 using U =
typename dolfinx::scalar_value_t<T>;
483 for (std::size_t i = 0; i < basis.size(); ++i)
487 V& bi = basis[i].get();
488 for (std::size_t j = 0; j < i; ++j)
490 const V& bj = basis[j].get();
494 std::ranges::transform(bj.array(), bi.array(), bi.array().begin(),
495 [dot_ij](
auto xj,
auto xi)
496 { return xi - dot_ij * xj; });
501 if (
norm *
norm < std::numeric_limits<U>::epsilon())
503 throw std::runtime_error(
504 "Linear dependency detected. Cannot orthogonalize.");
506 std::ranges::transform(bi.array(), bi.array().begin(),
507 [
norm](
auto x) { return x / norm; });
521 std::vector<std::reference_wrapper<const V>> basis,
522 dolfinx::scalar_value_t<typename V::value_type> eps = std::numeric_limits<
523 dolfinx::scalar_value_t<typename V::value_type>>::epsilon())
525 using T =
typename V::value_type;
526 for (std::size_t i = 0; i < basis.size(); i++)
528 for (std::size_t j = i; j < basis.size(); ++j)
530 T delta_ij = (i == j) ? T(1) : T(0);
532 if (std::norm(delta_ij - dot_ij) > eps)
container_type & mutable_array()
Get local part of the vector.
Definition Vector.h:353
container_type & array()
Get the process-local part of the vector.
Definition Vector.h:344
container_type::value_type value_type
Scalar type.
Definition Vector.h:118
const container_type & array() const
Get the process-local part of the vector (const version).
Definition Vector.h:349
void scatter_rev_begin(U pack, GetPtr get_ptr)
Start scatter (send) of ghost entry data to the owning process of an index.
Definition Vector.h:271
constexpr int bs() const
Get block size.
Definition Vector.h:339
Vector(std::shared_ptr< const common::IndexMap > map, int bs)
Create a distributed vector.
Definition Vector.h:128
void scatter_rev(BinaryOperation op)
Scatter (send) of ghost data values to the owning process and assign/accumulate into the owned data e...
Definition Vector.h:329
void scatter_fwd_end(U unpack)
End scatter (send) of local data values that are ghosted on other processes.
Definition Vector.h:213
void scatter_rev_end(U unpack)
End scatter of ghost data to owner and update owned entries.
Definition Vector.h:306
Vector & operator=(Vector &&x)=default
Move assignment operator.
std::shared_ptr< const common::IndexMap > index_map() const
Get IndexMap.
Definition Vector.h:336
void scatter_fwd_end()
End scatter (send) of local data values that are ghosted on other processes (simplified CPU version).
Definition Vector.h:229
void scatter_fwd()
Scatter (send) of local data values that are ghosted on other processes and update ghost entry values...
Definition Vector.h:246
Container container_type
Container type.
Definition Vector.h:115
void scatter_fwd_begin(U pack, GetPtr get_ptr)
Begin scatter (send) of local data that is ghosted on other processes.
Definition Vector.h:176
void set(value_type v)
Set all entries (including ghosts).
Definition Vector.h:154
Vector(Vector &&x)=default
Move constructor.
void scatter_fwd_begin()
Begin scatter (send) of local data that is ghosted on other processes (simplified CPU version).
Definition Vector.h:192
Vector(const Vector &x)=default
Copy constructor.
void scatter_rev_begin()
Start scatter (send) of ghost entry data to the owning process of an index (simplified CPU version).
Definition Vector.h:289
Access to pointer function concept.
Definition Vector.h:33
la::Vector scatter pack/unpack function concept.
Definition Vector.h:27
MPI_Datatype mpi_t
Retrieves the MPI data type associated to the provided type.
Definition MPI.h:280
Miscellaneous classes, functions and types.
Definition dolfinx_common.h:8
Linear algebra interface.
Definition sparsitybuild.h:15
void orthonormalize(std::vector< std::reference_wrapper< V > > basis)
Orthonormalize a set of vectors.
Definition Vector.h:477
auto squared_norm(const V &a)
Compute the squared L2 norm of vector.
Definition Vector.h:420
auto norm(const V &x, Norm type=Norm::l2)
Compute the norm of the vector.
Definition Vector.h:434
auto inner_product(const V &a, const V &b)
Compute the inner product of two vectors.
Definition Vector.h:389
bool is_orthonormal(std::vector< std::reference_wrapper< const V > > basis, dolfinx::scalar_value_t< typename V::value_type > eps=std::numeric_limits< dolfinx::scalar_value_t< typename V::value_type > >::epsilon())
Test if basis is orthonormal.
Definition Vector.h:520
Norm
Norm types.
Definition utils.h:17