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>);
 
   53  template <
typename, 
typename, 
typename>
 
   63    return [](
typename ScatterContainer::const_iterator idx_first,
 
   64              typename ScatterContainer::const_iterator idx_last,
 
   65              const auto in_first, 
auto out_first)
 
   68      std::transform(idx_first, idx_last, out_first,
 
   69                     [in_first](
auto p) { 
return *std::next(in_first, p); });
 
   79    return [](
typename ScatterContainer::const_iterator idx_first,
 
   80              typename ScatterContainer::const_iterator idx_last,
 
   81              const auto in_first, 
auto out_first)
 
   84      for (
typename ScatterContainer::const_iterator idx = idx_first;
 
   85           idx != idx_last; ++idx)
 
   87        std::size_t d = std::distance(idx_first, idx);
 
   88        *std::next(out_first, *idx) = *std::next(in_first, d);
 
   98  template <
typename BinaryOp>
 
   99  auto get_unpack_op(BinaryOp op)
 
  101    return [op](
typename ScatterContainer::const_iterator idx_first,
 
  102                typename ScatterContainer::const_iterator idx_last,
 
  103                const auto in_first, 
auto out_first)
 
  106      for (
typename ScatterContainer::const_iterator idx = idx_first;
 
  107           idx != idx_last; ++idx)
 
  109        std::size_t d = std::distance(idx_first, idx);
 
  110        auto& out = *std::next(out_first, *idx);
 
  111        out = op(out, *std::next(in_first, d));
 
  123  static_assert(std::is_same_v<value_type, typename container_type::value_type>,
 
  124                "Scalar type and container value type must be the same.");
 
  131  Vector(std::shared_ptr<const common::IndexMap> map, 
int bs)
 
  132      : _map(map), _bs(
bs), _x(
bs * (map->size_local() + map->num_ghosts())),
 
  134            std::make_shared<
common::Scatterer<ScatterContainer>>(*_map, 
bs)),
 
  135        _buffer_local(_scatterer->local_indices().size()),
 
  136        _buffer_remote(_scatterer->remote_indices().size())
 
 
  156  std::shared_ptr<const common::Scatterer<ScatterContainer>>
 
  157  scatter_ptr(
auto sc)
 const 
  159    using SC = 
typename std::remove_cv<
 
  160        typename decltype(sc)::element_type>::type::container_type;
 
  161    if constexpr (std::is_same_v<ScatterContainer, SC>)
 
  164      return std::make_shared<common::Scatterer<ScatterContainer>>(*sc);
 
  177  template <
typename T0, 
typename Container0, 
typename ScatterContainer0>
 
  178  explicit Vector(
const Vector<T0, Container0, ScatterContainer0>& x)
 
  179      : _map(x.
index_map()), _bs(x.
bs()), _x(x._x.begin(), x._x.end()),
 
  180        _scatterer(scatter_ptr(x._scatterer)), _request(MPI_REQUEST_NULL),
 
  181        _buffer_local(_scatterer->local_indices().size()),
 
  182        _buffer_remote(_scatterer->remote_indices().size())
 
 
  196  [[deprecated(
"Use std::ranges::fill(u.array(), v) instead.")]] 
void 
  199    std::ranges::fill(_x, v);
 
 
  216  template <
typename U, 
typename GetPtr>
 
  221    pack(_scatterer->local_indices().begin(), _scatterer->local_indices().end(),
 
  222         _x.begin(), _buffer_local.begin());
 
  223    _scatterer->scatter_fwd_begin(get_ptr(_buffer_local),
 
  224                                  get_ptr(_buffer_remote), _request);
 
 
  236    requires requires(Container c) {
 
  237      { c.data() } -> std::same_as<T*>;
 
 
  254  template <
typename U>
 
  255    requires VectorPackKernel<U, container_type, ScatterContainer>
 
  258    _scatterer->scatter_end(_request);
 
  259    unpack(_scatterer->remote_indices().begin(),
 
  260           _scatterer->remote_indices().end(), _buffer_remote.begin(),
 
  261           std::next(_x.begin(), _bs * _map->size_local()));
 
 
  273    requires requires(Container c) {
 
  274      { c.data() } -> std::same_as<T*>;
 
 
  290    requires requires(Container c) {
 
  291      { c.data() } -> std::same_as<T*>;
 
 
  311  template <
typename U, 
typename GetPtr>
 
  312    requires VectorPackKernel<U, container_type, ScatterContainer>
 
  313             and GetPtrConcept<GetPtr, Container, T>
 
  316    std::int32_t local_size = _bs * _map->size_local();
 
  317    pack(_scatterer->remote_indices().begin(),
 
  318         _scatterer->remote_indices().end(), std::next(_x.begin(), local_size),
 
  319         _buffer_remote.begin());
 
  320    _scatterer->scatter_rev_begin(get_ptr(_buffer_remote),
 
  321                                  get_ptr(_buffer_local), _request);
 
 
  333    requires requires(Container c) {
 
  334      { c.data() } -> std::same_as<T*>;
 
 
  347  template <
typename U>
 
  348    requires VectorPackKernel<U, container_type, ScatterContainer>
 
  351    _scatterer->scatter_end(_request);
 
  352    unpack(_scatterer->local_indices().begin(),
 
  353           _scatterer->local_indices().end(), _buffer_local.begin(),
 
 
  368  template <
class BinaryOperation>
 
  369    requires requires(Container c) {
 
  370      { c.data() } -> std::same_as<T*>;
 
  379  std::shared_ptr<const common::IndexMap> 
index_map()
 const { 
return _map; }
 
  382  constexpr int bs()
 const { 
return _bs; }
 
  403  std::shared_ptr<const common::IndexMap> _map;
 
  412  std::shared_ptr<const common::Scatterer<ScatterContainer>> _scatterer;
 
  415  MPI_Request _request = MPI_REQUEST_NULL;
 
 
  434  using T = 
typename V::value_type;
 
  435  const std::int32_t local_size = a.bs() * a.index_map()->size_local();
 
  436  if (local_size != b.bs() * b.index_map()->size_local())
 
  437    throw std::runtime_error(
"Incompatible vector sizes");
 
  439  const T local = std::transform_reduce(
 
  440      a.array().begin(), std::next(a.array().begin(), local_size),
 
  441      b.array().begin(), 
static_cast<T
>(0), std::plus{},
 
  444        if constexpr (std::is_same<T, std::complex<double>>::value
 
  445                      or std::is_same<T, std::complex<float>>::value)
 
  447          return std::conj(a) * b;
 
  455                a.index_map()->comm());
 
 
  465  using T = 
typename V::value_type;
 
  467  return std::real(result);
 
 
  479  using T = 
typename V::value_type;
 
  484    std::int32_t size_local = x.bs() * x.index_map()->size_local();
 
  485    using U = 
typename dolfinx::scalar_value_t<T>;
 
  486    U local_l1 = std::accumulate(
 
  487        x.array().begin(), std::next(x.array().begin(), size_local), U(0),
 
  488        [](
auto norm, 
auto x) { 
return norm + std::abs(x); });
 
  491                  x.index_map()->comm());
 
  498    std::int32_t size_local = x.bs() * x.index_map()->size_local();
 
  499    auto max_pos = std::max_element(
 
  500        x.array().begin(), std::next(x.array().begin(), size_local),
 
  501        [](T a, T b) { return std::norm(a) < std::norm(b); });
 
  502    auto local_linf = std::abs(*max_pos);
 
  503    decltype(local_linf) linf = 0;
 
  504    MPI_Allreduce(&local_linf, &linf, 1, 
MPI::mpi_t<
decltype(linf)>, MPI_MAX,
 
  505                  x.index_map()->comm());
 
  509    throw std::runtime_error(
"Norm type not supported");
 
 
  522  using T = 
typename V::value_type;
 
  523  using U = 
typename dolfinx::scalar_value_t<T>;
 
  526  for (std::size_t i = 0; i < basis.size(); ++i)
 
  530    V& bi = basis[i].get();
 
  531    for (std::size_t j = 0; j < i; ++j)
 
  533      const V& bj = basis[j].get();
 
  537      std::ranges::transform(bj.array(), bi.array(), bi.array().begin(),
 
  538                             [dot_ij](
auto xj, 
auto xi)
 
  539                             { return xi - dot_ij * xj; });
 
  544    if (
norm * 
norm < std::numeric_limits<U>::epsilon())
 
  546      throw std::runtime_error(
 
  547          "Linear dependency detected. Cannot orthogonalize.");
 
  549    std::ranges::transform(bi.array(), bi.array().begin(),
 
  550                           [
norm](
auto x) { return x / norm; });
 
 
  564    std::vector<std::reference_wrapper<const V>> basis,
 
  565    dolfinx::scalar_value_t<typename V::value_type> eps = std::numeric_limits<
 
  566        dolfinx::scalar_value_t<typename V::value_type>>::epsilon())
 
  568  using T = 
typename V::value_type;
 
  569  for (std::size_t i = 0; i < basis.size(); i++)
 
  571    for (std::size_t j = i; j < basis.size(); ++j)
 
  573      T delta_ij = (i == j) ? T(1) : T(0);
 
  575      if (std::norm(delta_ij - dot_ij) > eps)
 
 
A vector that can be distributed across processes.
Definition Vector.h:50
 
container_type & mutable_array()
Get local part of the vector.
Definition Vector.h:396
 
container_type & array()
Get the process-local part of the vector.
Definition Vector.h:387
 
container_type::value_type value_type
Scalar type.
Definition Vector.h:121
 
const container_type & array() const
Get the process-local part of the vector (const version).
Definition Vector.h:392
 
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:314
 
Vector(const Vector< T0, Container0, ScatterContainer0 > &x)
Copy-convert vector, possibly using to different container types.
Definition Vector.h:178
 
constexpr int bs() const
Get block size.
Definition Vector.h:382
 
Vector(std::shared_ptr< const common::IndexMap > map, int bs)
Create a distributed vector.
Definition Vector.h:131
 
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:372
 
void scatter_fwd_end(U unpack)
End scatter (send) of local data values that are ghosted on other processes.
Definition Vector.h:256
 
void scatter_rev_end(U unpack)
End scatter of ghost data to owner and update owned entries.
Definition Vector.h:349
 
Vector & operator=(Vector &&x)=default
Move assignment operator.
 
std::shared_ptr< const common::IndexMap > index_map() const
Get IndexMap.
Definition Vector.h:379
 
void scatter_fwd_end()
End scatter (send) of local data values that are ghosted on other processes (simplified CPU version).
Definition Vector.h:272
 
void scatter_fwd()
Scatter (send) of local data values that are ghosted on other processes and update ghost entry values...
Definition Vector.h:289
 
Container container_type
Container type.
Definition Vector.h:118
 
void scatter_fwd_begin(U pack, GetPtr get_ptr)
Begin scatter (send) of local data that is ghosted on other processes.
Definition Vector.h:219
 
void set(value_type v)
Set all entries (including ghosts).
Definition Vector.h:197
 
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:235
 
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:332
 
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:520
 
auto squared_norm(const V &a)
Compute the squared L2 norm of vector.
Definition Vector.h:463
 
auto norm(const V &x, Norm type=Norm::l2)
Compute the norm of the vector.
Definition Vector.h:477
 
auto inner_product(const V &a, const V &b)
Compute the inner product of two vectors.
Definition Vector.h:432
 
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:563
 
Norm
Norm types.
Definition utils.h:17