libMesh
distributed_vector.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // Local includes
21 #include "libmesh/distributed_vector.h"
22 
23 // libMesh includes
24 #include "libmesh/dense_vector.h"
25 #include "libmesh/dense_subvector.h"
26 #include "libmesh/int_range.h"
27 #include "libmesh/libmesh_common.h"
28 #include "libmesh/tensor_tools.h"
29 
30 // TIMPI includes
32 #include "timpi/parallel_sync.h"
33 
34 // C++ includes
35 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
36 #include <cmath> // for std::abs
37 #include <limits> // std::numeric_limits<T>::min()
38 
39 
40 namespace libMesh
41 {
42 
43 
44 
45 //--------------------------------------------------------------------------
46 // DistributedVector methods
47 template <typename T>
49 {
50  // This function must be run on all processors at once
51  parallel_object_only();
52 
53  libmesh_assert (this->initialized());
54  libmesh_assert_equal_to (_values.size(), _local_size);
55  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
56 
57  T local_sum = 0.;
58 
59  for (auto & val : _values)
60  local_sum += val;
61 
62  this->comm().sum(local_sum);
63 
64  return local_sum;
65 }
66 
67 
68 
69 template <typename T>
71 {
72  // This function must be run on all processors at once
73  parallel_object_only();
74 
75  libmesh_assert (this->initialized());
76  libmesh_assert_equal_to (_values.size(), _local_size);
77  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
78 
79  Real local_l1 = 0.;
80 
81  for (auto & val : _values)
82  local_l1 += std::abs(val);
83 
84  this->comm().sum(local_l1);
85 
86  return local_l1;
87 }
88 
89 
90 
91 template <typename T>
93 {
94  // This function must be run on all processors at once
95  parallel_object_only();
96 
97  libmesh_assert (this->initialized());
98  libmesh_assert_equal_to (_values.size(), _local_size);
99  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
100 
101  Real local_l2 = 0.;
102 
103  for (auto & val : _values)
104  local_l2 += TensorTools::norm_sq(val);
105 
106  this->comm().sum(local_l2);
107 
108  return std::sqrt(local_l2);
109 }
110 
111 
112 
113 template <typename T>
115 {
116  // This function must be run on all processors at once
117  parallel_object_only();
118 
119  libmesh_assert (this->initialized());
120  libmesh_assert_equal_to (_values.size(), _local_size);
121  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
122 
123  Real local_linfty = 0.;
124 
125  for (auto & val : _values)
126  local_linfty = std::max(local_linfty,
127  static_cast<Real>(std::abs(val))
128  ); // Note we static_cast so that both
129  // types are the same, as required
130  // by std::max
131 
132  this->comm().max(local_linfty);
133 
134  return local_linfty;
135 }
136 
137 
138 
139 template <typename T>
141 {
142  libmesh_assert (this->closed());
143  libmesh_assert (this->initialized());
144  libmesh_assert_equal_to (_values.size(), _local_size);
145  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
146 
147  add(1., v);
148 
149  return *this;
150 }
151 
152 
153 
154 template <typename T>
156 {
157  libmesh_assert (this->closed());
158  libmesh_assert (this->initialized());
159  libmesh_assert_equal_to (_values.size(), _local_size);
160  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
161 
162  add(-1., v);
163 
164  return *this;
165 }
166 
167 
168 
169 template <typename T>
171 {
172  libmesh_assert_equal_to(size(), v.size());
173 
174  const DistributedVector<T> & v_vec = cast_ref<const DistributedVector<T> &>(v);
175 
176  for (auto i : index_range(_values))
177  _values[i] *= v_vec._values[i];
178 
179  return *this;
180 }
181 
182 
183 
184 template <typename T>
186 {
187  libmesh_assert_equal_to(size(), v.size());
188 
189  const DistributedVector<T> & v_vec = cast_ref<const DistributedVector<T> &>(v);
190 
191  for (auto i : index_range(_values))
192  _values[i] /= v_vec._values[i];
193 
194  return *this;
195 }
196 
197 
198 
199 
200 template <typename T>
202 {
203  for (auto & val : _values)
204  {
205  // Don't divide by zero
206  libmesh_assert_not_equal_to (val, T(0));
207 
208  val = 1. / val;
209  }
210 }
211 
212 
213 
214 
215 template <typename T>
217 {
218  // Replace values by complex conjugate
219  for (auto & val : _values)
220  val = libmesh_conj(val);
221 }
222 
223 
224 
225 
226 
227 template <typename T>
229 {
230  libmesh_assert (this->initialized());
231  libmesh_assert_equal_to (_values.size(), _local_size);
232  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
233 
234  for (auto & val : _values)
235  val += v;
236 }
237 
238 
239 
240 template <typename T>
242 {
243  libmesh_assert (this->initialized());
244  libmesh_assert_equal_to (_values.size(), _local_size);
245  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
246 
247  add (1., v);
248 }
249 
250 
251 
252 template <typename T>
253 void DistributedVector<T>::add (const T a, const NumericVector<T> & v_in)
254 {
255  libmesh_assert (this->initialized());
256  libmesh_assert_equal_to (_values.size(), _local_size);
257  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
258 
259  // Make sure the NumericVector passed in is really a DistributedVector
260  const DistributedVector<T> * v = cast_ptr<const DistributedVector<T> *>(&v_in);
261  if (!v)
262  libmesh_error_msg("Cannot add different types of NumericVectors.");
263 
264  for (auto i : index_range(_values))
265  _values[i] += a * v->_values[i];
266 }
267 
268 
269 
270 template <typename T>
271 void DistributedVector<T>::scale (const T factor)
272 {
273  libmesh_assert (this->initialized());
274  libmesh_assert_equal_to (_values.size(), _local_size);
275  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
276 
277  for (auto & val : _values)
278  val *= factor;
279 }
280 
281 template <typename T>
283 {
284  libmesh_assert (this->initialized());
285  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
286 
287  for (auto & val : _values)
288  val = std::abs(val);
289 }
290 
291 
292 
293 
294 
295 template <typename T>
297 {
298  // This function must be run on all processors at once
299  parallel_object_only();
300 
301  // Make sure the NumericVector passed in is really a DistributedVector
302  const DistributedVector<T> * v = cast_ptr<const DistributedVector<T> *>(&V);
303 
304  // Make sure that the two vectors are distributed in the same way.
305  libmesh_assert_equal_to ( this->first_local_index(), v->first_local_index() );
306  libmesh_assert_equal_to ( this->last_local_index(), v->last_local_index() );
307 
308  // The result of dotting together the local parts of the vector.
309  T local_dot = 0;
310 
311  for (auto i : index_range(_values))
312  local_dot += this->_values[i] * v->_values[i];
313 
314  // The local dot products are now summed via MPI
315  this->comm().sum(local_dot);
316 
317  return local_dot;
318 }
319 
320 
321 
322 template <typename T>
325 {
326  libmesh_assert (this->initialized());
327  libmesh_assert_equal_to (_values.size(), _local_size);
328  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
329 
330  for (auto & val : _values)
331  val = s;
332 
333  return *this;
334 }
335 
336 
337 
338 template <typename T>
341 {
342  // Make sure the NumericVector passed in is really a DistributedVector
343  const DistributedVector<T> * v = cast_ptr<const DistributedVector<T> *>(&v_in);
344 
345  *this = *v;
346 
347  return *this;
348 }
349 
350 
351 
352 template <typename T>
355 {
357  this->_is_closed = v._is_closed;
358 
359  _global_size = v._global_size;
360  _local_size = v._local_size;
361  _first_local_index = v._first_local_index;
362  _last_local_index = v._last_local_index;
363 
364  if (v.local_size() == this->local_size())
365  _values = v._values;
366 
367  else
368  libmesh_error_msg("v.local_size() = " << v.local_size() << " must be equal to this->local_size() = " << this->local_size());
369 
370  return *this;
371 }
372 
373 
374 
375 template <typename T>
377 DistributedVector<T>::operator = (const std::vector<T> & v)
378 {
379  libmesh_assert (this->initialized());
380  libmesh_assert_equal_to (_values.size(), _local_size);
381  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
382 
383  if (v.size() == local_size())
384  _values = v;
385 
386  else if (v.size() == size())
387  for (auto i : index_range(*this))
388  _values[i-first_local_index()] = v[i];
389 
390  else
391  libmesh_error_msg("Incompatible sizes in DistributedVector::operator=");
392 
393  return *this;
394 }
395 
396 
397 
398 template <typename T>
400 
401 {
402  libmesh_assert (this->initialized());
403  libmesh_assert_equal_to (_values.size(), _local_size);
404  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
405 
406  DistributedVector<T> * v_local = cast_ptr<DistributedVector<T> *>(&v_local_in);
407 
408  v_local->_first_local_index = 0;
409 
410  v_local->_global_size =
411  v_local->_local_size =
412  v_local->_last_local_index = size();
413 
414  v_local->_is_initialized =
415  v_local->_is_closed = true;
416 
417  // Call localize on the vector's values. This will help
418  // prevent code duplication
419  localize (v_local->_values);
420 
421 #ifndef LIBMESH_HAVE_MPI
422 
423  libmesh_assert_equal_to (local_size(), size());
424 
425 #endif
426 }
427 
428 
429 
430 template <typename T>
432  const std::vector<numeric_index_type> &) const
433 {
434  libmesh_assert (this->initialized());
435  libmesh_assert_equal_to (_values.size(), _local_size);
436  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
437 
438  // TODO: We don't yet support the send list; this is inefficient:
439  localize (v_local_in);
440 }
441 
442 
443 
444 template <typename T>
445 void DistributedVector<T>::localize (std::vector<T> & v_local,
446  const std::vector<numeric_index_type> & indices) const
447 {
448  // Resize v_local so there is enough room to hold all the local values.
449  v_local.resize(indices.size());
450 
451  // We need to know who has the values we want, so get everyone's _local_size
452  std::vector<numeric_index_type> local_sizes;
453  this->comm().allgather (_local_size, local_sizes);
454 
455  // Make a vector of partial sums of local sizes
456  std::vector<numeric_index_type> local_size_sums(this->n_processors());
457  local_size_sums[0] = local_sizes[0];
458  for (auto i : IntRange<numeric_index_type>(1, local_sizes.size()))
459  local_size_sums[i] = local_size_sums[i-1] + local_sizes[i];
460 
461  // We now fill in 'requested_ids' based on the indices. Also keep
462  // track of the local index (in the indices vector) for each of
463  // these, since we need that when unpacking.
464  std::map<processor_id_type, std::vector<numeric_index_type>>
465  requested_ids, local_requested_ids;
466 
467  // We'll use this typedef a couple of times below.
468  typedef typename std::vector<numeric_index_type>::iterator iter_t;
469 
470  // For each index in indices, determine which processor it is on.
471  // This is an O(N*log(p)) algorithm that uses std::upper_bound().
472  // Note: upper_bound() returns an iterator to the first entry which is
473  // greater than the given value.
474  for (auto i : index_range(indices))
475  {
476  iter_t ub = std::upper_bound(local_size_sums.begin(),
477  local_size_sums.end(),
478  indices[i]);
479 
480  processor_id_type on_proc = cast_int<processor_id_type>
481  (std::distance(local_size_sums.begin(), ub));
482 
483  requested_ids[on_proc].push_back(indices[i]);
484  local_requested_ids[on_proc].push_back(i);
485  }
486 
487  auto gather_functor =
488  [this]
489  (processor_id_type, const std::vector<dof_id_type> & ids,
490  std::vector<T> & values)
491  {
492  // The first send/receive we did was for indices, the second one will be
493  // for corresponding floating point values, so create storage for that now...
494  const std::size_t ids_size = ids.size();
495  values.resize(ids_size);
496 
497  for (std::size_t i=0; i != ids_size; i++)
498  {
499  // The index of the requested value
500  const numeric_index_type requested_index = ids[i];
501 
502  // Transform into local numbering, and get requested value.
503  values[i] = _values[requested_index - _first_local_index];
504  }
505  };
506 
507  auto action_functor =
508  [& v_local, & local_requested_ids]
509  (processor_id_type pid,
510  const std::vector<dof_id_type> &,
511  const std::vector<T> & values)
512  {
513  // Now write the received values to the appropriate place(s) in v_local
514  for (auto i : index_range(values))
515  {
516  libmesh_assert(local_requested_ids.count(pid));
517  libmesh_assert_less(i, local_requested_ids[pid].size());
518 
519  // Get the index in v_local where this value needs to be inserted.
520  const numeric_index_type local_requested_index =
521  local_requested_ids[pid][i];
522 
523  // Actually set the value in v_local
524  v_local[local_requested_index] = values[i];
525  }
526  };
527 
528  const T * ex = nullptr;
529  Parallel::pull_parallel_vector_data
530  (this->comm(), requested_ids, gather_functor, action_functor, ex);
531 }
532 
533 
534 
535 template <typename T>
537  const numeric_index_type last_local_idx,
538  const std::vector<numeric_index_type> & send_list)
539 {
540  // Only good for serial vectors
541  libmesh_assert_equal_to (this->size(), this->local_size());
542  libmesh_assert_greater (last_local_idx, first_local_idx);
543  libmesh_assert_less_equal (send_list.size(), this->size());
544  libmesh_assert_less (last_local_idx, this->size());
545 
546  const numeric_index_type my_size = this->size();
547  const numeric_index_type my_local_size = (last_local_idx - first_local_idx + 1);
548 
549  // Don't bother for serial cases
550  if ((first_local_idx == 0) &&
551  (my_local_size == my_size))
552  return;
553 
554 
555  // Build a parallel vector, initialize it with the local
556  // parts of (*this)
557  DistributedVector<T> parallel_vec(this->comm());
558 
559  parallel_vec.init (my_size, my_local_size, true, PARALLEL);
560 
561  // Copy part of *this into the parallel_vec
562  for (numeric_index_type i=first_local_idx; i<=last_local_idx; i++)
563  parallel_vec._values[i-first_local_idx] = _values[i];
564 
565  // localize like normal
566  parallel_vec.localize (*this, send_list);
567 }
568 
569 
570 
571 template <typename T>
572 void DistributedVector<T>::localize (std::vector<T> & v_local) const
573 {
574  // This function must be run on all processors at once
575  parallel_object_only();
576 
577  libmesh_assert (this->initialized());
578  libmesh_assert_equal_to (_values.size(), _local_size);
579  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
580 
581  v_local = this->_values;
582 
583  this->comm().allgather (v_local);
584 
585 #ifndef LIBMESH_HAVE_MPI
586  libmesh_assert_equal_to (local_size(), size());
587 #endif
588 }
589 
590 
591 
592 template <typename T>
593 void DistributedVector<T>::localize_to_one (std::vector<T> & v_local,
594  const processor_id_type pid) const
595 {
596  // This function must be run on all processors at once
597  parallel_object_only();
598 
599  libmesh_assert (this->initialized());
600  libmesh_assert_equal_to (_values.size(), _local_size);
601  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
602 
603  v_local = this->_values;
604 
605  this->comm().gather (pid, v_local);
606 
607 #ifndef LIBMESH_HAVE_MPI
608  libmesh_assert_equal_to (local_size(), size());
609 #endif
610 }
611 
612 
613 
614 template <typename T>
616  const NumericVector<T> &)
617 //void DistributedVector<T>::pointwise_mult (const NumericVector<T> & vec1,
618 // const NumericVector<T> & vec2)
619 {
620  libmesh_not_implemented();
621 }
622 
623 
624 
625 //--------------------------------------------------------------
626 // Explicit instantiations
627 template class DistributedVector<Number>;
628 
629 } // namespace libMesh
libMesh::PARALLEL
Definition: enum_parallel_type.h:36
libMesh::DistributedVector::operator*=
virtual NumericVector< T > & operator*=(const NumericVector< T > &v) override
Computes the component-wise multiplication of this vector's entries by another's, .
Definition: distributed_vector.C:170
libMesh::DistributedVector::conjugate
virtual void conjugate() override
Negates the imaginary component of each entry in the vector.
Definition: distributed_vector.C:216
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
std::sqrt
MetaPhysicL::DualNumber< T, D > sqrt(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::DistributedVector::pointwise_mult
virtual void pointwise_mult(const NumericVector< T > &vec1, const NumericVector< T > &vec2) override
Computes (summation not implied) i.e.
Definition: distributed_vector.C:615
libMesh::libmesh_conj
T libmesh_conj(T a)
Definition: libmesh_common.h:167
libMesh::NumericVector::_is_initialized
bool _is_initialized
true once init() has been called.
Definition: numeric_vector.h:732
libMesh::NumericVector::_is_closed
bool _is_closed
Flag which tracks whether the vector's values are consistent on all processors after insertion or add...
Definition: numeric_vector.h:727
libMesh::DistributedVector::localize
virtual void localize(std::vector< T > &v_local) const override
Creates a copy of the global vector in the local vector v_local.
Definition: distributed_vector.C:572
libMesh::DistributedVector::_local_size
numeric_index_type _local_size
The local vector size.
Definition: distributed_vector.h:245
libMesh::DistributedVector::_values
std::vector< T > _values
Actual vector datatype to hold vector entries.
Definition: distributed_vector.h:235
libMesh::DistributedVector::localize_to_one
virtual void localize_to_one(std::vector< T > &v_local, const processor_id_type proc_id=0) const override
Creates a local copy of the global vector in v_local only on processor proc_id.
Definition: distributed_vector.C:593
libMesh::DistributedVector::first_local_index
virtual numeric_index_type first_local_index() const override
Definition: distributed_vector.h:526
parallel_sync.h
libMesh::DistributedVector::init
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) override
Change the dimension of the vector to n.
Definition: distributed_vector.h:317
libMesh::DistributedVector::l2_norm
virtual Real l2_norm() const override
Definition: distributed_vector.C:92
libMesh::DistributedVector::abs
virtual void abs() override
Sets for each entry in the vector.
Definition: distributed_vector.C:282
libMesh::NumericVector::size
virtual numeric_index_type size() const =0
libMesh::NumericVector
Provides a uniform interface to vector storage schemes for different linear algebra libraries.
Definition: vector_fe_ex5.C:43
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::IntRange
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
Definition: int_range.h:53
libMesh::DistributedVector::sum
virtual T sum() const override
Definition: distributed_vector.C:48
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
parallel_implementation.h
libMesh::DistributedVector::_first_local_index
numeric_index_type _first_local_index
The first component stored locally.
Definition: distributed_vector.h:250
libMesh::processor_id_type
uint8_t processor_id_type
Definition: id_types.h:104
libMesh::DistributedVector::last_local_index
virtual numeric_index_type last_local_index() const override
Definition: distributed_vector.h:539
libMesh::numeric_index_type
dof_id_type numeric_index_type
Definition: id_types.h:99
libMesh::initialized
bool initialized()
Checks that library initialization has been done.
Definition: libmesh.C:265
libMesh::DistributedVector::_last_local_index
numeric_index_type _last_local_index
The last component (+1) stored locally.
Definition: distributed_vector.h:255
distance
Real distance(const Point &p)
Definition: subdomains_ex3.C:50
libMesh::DistributedVector::dot
virtual T dot(const NumericVector< T > &V) const override
Definition: distributed_vector.C:296
libMesh::DistributedVector::scale
virtual void scale(const T factor) override
Scale each element of the vector by the given factor.
Definition: distributed_vector.C:271
libMesh::DistributedVector
This class provides a simple parallel, distributed vector datatype which is specific to libmesh.
Definition: distributed_vector.h:53
libMesh::TensorTools::norm_sq
T norm_sq(std::complex< T > a)
Definition: tensor_tools.h:85
libMesh::DistributedVector::operator/=
virtual NumericVector< T > & operator/=(const NumericVector< T > &v) override
Computes the component-wise division of this vector's entries by another's, .
Definition: distributed_vector.C:185
libMesh::libMeshPrivateData::_is_initialized
bool _is_initialized
Flag that tells if init() has been called.
Definition: libmesh.C:246
libMesh::DistributedVector::operator=
DistributedVector & operator=(const DistributedVector &)
Copy assignment operator.
Definition: distributed_vector.C:354
libMesh::DistributedVector::add
virtual void add(const numeric_index_type i, const T value) override
Adds value to each entry of the vector.
Definition: distributed_vector.h:582
libMesh::DistributedVector::linfty_norm
virtual Real linfty_norm() const override
Definition: distributed_vector.C:114
libMesh::DistributedVector::_global_size
numeric_index_type _global_size
The global vector size.
Definition: distributed_vector.h:240
libMesh::DistributedVector::operator-=
virtual NumericVector< T > & operator-=(const NumericVector< T > &v) override
Subtracts v from *this, .
Definition: distributed_vector.C:155
libMesh::DistributedVector::local_size
virtual numeric_index_type local_size() const override
Definition: distributed_vector.h:513
libMesh::DistributedVector::operator+=
virtual NumericVector< T > & operator+=(const NumericVector< T > &v) override
Adds v to *this, .
Definition: distributed_vector.C:140
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::DistributedVector::reciprocal
virtual void reciprocal() override
Computes the component-wise reciprocal, .
Definition: distributed_vector.C:201
libMesh::DistributedVector::l1_norm
virtual Real l1_norm() const override
Definition: distributed_vector.C:70
libMesh::closed
bool closed()
Checks that the library has been closed.
Definition: libmesh.C:272