12#include <dolfinx/common/IndexMap.h> 
   13#include <dolfinx/common/Scatterer.h> 
   14#include <dolfinx/common/types.h> 
   29template <
typename T, 
typename Container = std::vector<T>>
 
   32  static_assert(std::is_same_v<typename Container::value_type, T>);
 
   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.");
 
   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()))
 
 
   57      : _map(x._map), _scatterer(x._scatterer), _bs(x._bs),
 
   59        _buffer_remote(x._buffer_remote), _x(x._x)
 
 
   65      : _map(std::
move(x._map)), _scatterer(std::
move(x._scatterer)),
 
   66        _bs(std::
move(x._bs)),
 
   68        _buffer_local(std::move(x._buffer_local)),
 
   69        _buffer_remote(std::move(x._buffer_remote)), _x(std::move(x._x))
 
 
   74  Vector& operator=(
const Vector& x) = 
delete;
 
   87    const std::int32_t local_size = _bs * _map->size_local();
 
   88    std::span<const value_type> 
x_local(_x.data(), local_size);
 
   90    auto pack = [](
const auto& 
in, 
const auto& 
idx, 
auto& 
out)
 
   92      for (std::size_t 
i = 0; 
i < 
idx.size(); ++
i)
 
   95    pack(
x_local, _scatterer->local_indices(), _buffer_local);
 
   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));
 
 
  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));
 
  113      for (std::size_t 
i = 0; 
i < 
idx.size(); ++
i)
 
  118           [](
auto , 
auto b) { return b; });
 
 
  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);
 
  137    auto pack = [](
const auto& 
in, 
const auto& 
idx, 
auto& 
out)
 
  139      for (std::size_t 
i = 0; 
i < 
idx.size(); ++
i)
 
  142    pack(
x_remote, _scatterer->remote_indices(), _buffer_remote);
 
  144    _scatterer->scatter_rev_begin(std::span<const value_type>(_buffer_remote),
 
  145                                  std::span<value_type>(_buffer_local),
 
 
  155  template <
class BinaryOperation>
 
  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);
 
  164      for (std::size_t 
i = 0; 
i < 
idx.size(); ++
i)
 
 
  175  template <
class BinaryOperation>
 
  183  std::shared_ptr<const common::IndexMap> 
index_map()
 const { 
return _map; }
 
  186  constexpr int bs()
 const { 
return _bs; }
 
  189  std::span<const value_type> 
array()
 const 
  191    return std::span<const value_type>(_x);
 
 
  199  std::shared_ptr<const common::IndexMap> _map;
 
  202  std::shared_ptr<const common::Scatterer<>> _scatterer;
 
 
  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);
 
  233  const T local = std::transform_reduce(
 
  234      x_a.begin(), x_a.end(), x_b.begin(), 
static_cast<T
>(0), std::plus{},
 
  237        if constexpr (std::is_same<T, std::complex<double>>::value
 
  238                      or std::is_same<T, std::complex<float>>::value)
 
  240          return std::conj(a) * b;
 
  247  MPI_Allreduce(&local, &result, 1, dolfinx::MPI::mpi_type<T>(), MPI_SUM,
 
  248                a.index_map()->comm());
 
 
  257  using T = 
typename V::value_type;
 
  259  return std::real(result);
 
 
  269  using T = 
typename V::value_type;
 
  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>;
 
  278        = std::accumulate(data.begin(), data.end(), U(0),
 
  279                          [](
auto norm, 
auto x) { 
return norm + std::abs(x); });
 
  281    MPI_Allreduce(&local_l1, &l1, 1, MPI::mpi_type<U>(), MPI_SUM,
 
  282                  x.index_map()->comm());
 
  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(),
 
  293                                    { return std::norm(a) < std::norm(b); });
 
  294    auto local_linf = std::abs(*max_pos);
 
  295    decltype(local_linf) linf = 0;
 
  296    MPI_Allreduce(&local_linf, &linf, 1, MPI::mpi_type<
decltype(linf)>(),
 
  297                  MPI_MAX, x.index_map()->comm());
 
  301    throw std::runtime_error(
"Norm type not supported");
 
 
  313  using T = 
typename V::value_type;
 
  314  using U = 
typename dolfinx::scalar_value_type_t<T>;
 
  317  for (std::size_t i = 0; i < basis.size(); ++i)
 
  321    V& bi = basis[i].get();
 
  322    for (std::size_t j = 0; j < i; ++j)
 
  324      const V& bj = basis[j].get();
 
  328      std::transform(bj.array().begin(), bj.array().end(), bi.array().begin(),
 
  329                     bi.mutable_array().begin(),
 
  330                     [dot_ij](
auto xj, 
auto xi) { return xi - dot_ij * xj; });
 
  334    auto norm = la::norm(bi, la::Norm::l2);
 
  335    if (
norm * 
norm < std::numeric_limits<U>::epsilon())
 
  337      throw std::runtime_error(
 
  338          "Linear dependency detected. Cannot orthogonalize.");
 
  340    std::transform(bi.array().begin(), bi.array().end(),
 
  341                   bi.mutable_array().begin(),
 
  342                   [
norm](
auto x) { return x / norm; });
 
 
  352  using T = 
typename V::value_type;
 
  353  using U = 
typename dolfinx::scalar_value_type_t<T>;
 
  356  for (std::size_t i = 0; i < basis.size(); i++)
 
  358    for (std::size_t j = i; j < basis.size(); j++)
 
  360      T delta_ij = (i == j) ? T(1) : T(0);
 
  362      if (std::norm(delta_ij - dot_ij) > std::numeric_limits<U>::epsilon())
 
 
Distributed vector.
Definition Vector.h:31
 
void scatter_fwd_begin()
Begin scatter of local data from owner to ghosts on other ranks.
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)
Create a distributed vector.
Definition Vector.h:47
 
Vector & operator=(Vector &&x)=default
Move Assignment operator.
 
void scatter_fwd()
Scatter local data to ghost positions on other ranks.
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()
End scatter of local data from owner to ghosts on other ranks.
Definition Vector.h:104
 
void scatter_rev_begin()
Start scatter of ghost data to owner.
Definition Vector.h:131
 
void set(value_type v)
Set all entries (including ghosts)
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)
End scatter of ghost data to owner. This process may receive data from more than one process,...
Definition Vector.h:156
 
void scatter_rev(BinaryOperation op)
Scatter ghost data to owner. This process may receive data from more than one process,...
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)
Orthonormalize a set of vectors.
Definition Vector.h:311
 
Norm
Norm types.
Definition utils.h:17
 
auto squared_norm(const V &a)
Compute the squared L2 norm of vector.
Definition Vector.h:255
 
bool is_orthonormal(std::vector< std::reference_wrapper< const V > > basis)
Test if basis is orthonormal.
Definition Vector.h:350
 
auto norm(const V &x, Norm type=Norm::l2)
Compute the norm of the vector.
Definition Vector.h:267
 
auto inner_product(const V &a, const V &b)
Compute the inner product of two vectors. The two vectors must have the same parallel layout.
Definition Vector.h:224