$extrastylesheet
distributed_vector.h
Go to the documentation of this file.
00001 // The libMesh Finite Element Library.
00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
00003 
00004 // This library is free software; you can redistribute it and/or
00005 // modify it under the terms of the GNU Lesser General Public
00006 // License as published by the Free Software Foundation; either
00007 // version 2.1 of the License, or (at your option) any later version.
00008 
00009 // This library is distributed in the hope that it will be useful,
00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00012 // Lesser General Public License for more details.
00013 
00014 // You should have received a copy of the GNU Lesser General Public
00015 // License along with this library; if not, write to the Free Software
00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00017 
00018 
00019 
00020 #include "libmesh/libmesh_common.h"
00021 
00022 
00023 
00024 #ifndef LIBMESH_DISTRIBUTED_VECTOR_H
00025 #define LIBMESH_DISTRIBUTED_VECTOR_H
00026 
00027 // Local includes
00028 #include "libmesh/numeric_vector.h"
00029 #include "libmesh/parallel.h"
00030 
00031 // C++ includes
00032 #include <vector>
00033 #include <algorithm>
00034 #include <limits>
00035 
00036 namespace libMesh
00037 {
00038 
00039 
00040 
00051 template <typename T>
00052 class DistributedVector : public NumericVector<T>
00053 {
00054 public:
00055 
00059   explicit
00060   DistributedVector (const Parallel::Communicator &comm,
00061                      const ParallelType = AUTOMATIC);
00062 
00066   explicit
00067   DistributedVector (const Parallel::Communicator &comm,
00068                      const numeric_index_type n,
00069                      const ParallelType ptype = AUTOMATIC);
00070 
00075   DistributedVector (const Parallel::Communicator &comm,
00076                      const numeric_index_type n,
00077                      const numeric_index_type n_local,
00078                      const ParallelType ptype = AUTOMATIC);
00079 
00085   DistributedVector (const Parallel::Communicator &comm,
00086                      const numeric_index_type N,
00087                      const numeric_index_type n_local,
00088                      const std::vector<numeric_index_type>& ghost,
00089                      const ParallelType ptype = AUTOMATIC);
00090 
00095   ~DistributedVector ();
00096 
00100   void close ();
00101 
00105   void clear ();
00106 
00111   void zero ();
00112 
00118   virtual UniquePtr<NumericVector<T> > zero_clone () const;
00119 
00123   UniquePtr<NumericVector<T> > clone () const;
00124 
00137   void init (const numeric_index_type N,
00138              const numeric_index_type n_local,
00139              const bool         fast=false,
00140              const ParallelType ptype=AUTOMATIC);
00141 
00145   void init (const numeric_index_type N,
00146              const bool         fast=false,
00147              const ParallelType ptype=AUTOMATIC);
00148 
00153   virtual void init (const numeric_index_type /*N*/,
00154                      const numeric_index_type /*n_local*/,
00155                      const std::vector<numeric_index_type>& /*ghost*/,
00156                      const bool /*fast*/ = false,
00157                      const ParallelType = AUTOMATIC);
00158 
00163   virtual void init (const NumericVector<T>& other,
00164                      const bool fast = false);
00165 
00169   NumericVector<T> & operator= (const T s);
00170 
00174   NumericVector<T> & operator= (const NumericVector<T> &V);
00175 
00179   DistributedVector<T> & operator= (const DistributedVector<T> &V);
00180 
00184   NumericVector<T> & operator= (const std::vector<T> &v);
00185 
00191   Real min () const;
00192 
00198   Real max () const;
00199 
00203   T sum() const;
00204 
00209   Real l1_norm () const;
00210 
00216   Real l2_norm () const;
00217 
00223   Real linfty_norm () const;
00224 
00232   numeric_index_type size () const;
00233 
00238   numeric_index_type local_size() const;
00239 
00244   numeric_index_type first_local_index() const;
00245 
00250   numeric_index_type last_local_index() const;
00251 
00255   T operator() (const numeric_index_type i) const;
00256 
00261   NumericVector<T> & operator += (const NumericVector<T> &V);
00262 
00267   NumericVector<T> & operator -= (const NumericVector<T> &V);
00268 
00272   virtual NumericVector<T> & operator /= (NumericVector<T> & v);
00273 
00277   virtual void reciprocal();
00278 
00283   virtual void conjugate();
00284 
00288   void set (const numeric_index_type i, const T value);
00289 
00293   void add (const numeric_index_type i, const T value);
00294 
00300   void add (const T s);
00301 
00307   void add (const NumericVector<T>& V);
00308 
00314   void add (const T a, const NumericVector<T>& v);
00315 
00320   using NumericVector<T>::add_vector;
00321 
00328   void add_vector (const NumericVector<T>&,
00329                    const SparseMatrix<T>&)
00330   { libmesh_not_implemented(); }
00331 
00338   void add_vector_transpose (const NumericVector<T>&,
00339                              const SparseMatrix<T>&)
00340   { libmesh_not_implemented(); }
00341 
00346   void scale (const T factor);
00347 
00352   virtual void abs();
00353 
00357   virtual T dot(const NumericVector<T>& V) const;
00358 
00363   void localize (std::vector<T>& v_local) const;
00364 
00369   void localize (NumericVector<T>& v_local) const;
00370 
00376   void localize (NumericVector<T>& v_local,
00377                  const std::vector<numeric_index_type>& send_list) const;
00378 
00383   void localize (const numeric_index_type first_local_idx,
00384                  const numeric_index_type last_local_idx,
00385                  const std::vector<numeric_index_type>& send_list);
00386 
00393   void localize_to_one (std::vector<T>& v_local,
00394                         const processor_id_type proc_id=0) const;
00395 
00400   virtual void pointwise_mult (const NumericVector<T>& vec1,
00401                                const NumericVector<T>& vec2);
00402 
00406   virtual void swap (NumericVector<T> &v);
00407 
00408 private:
00409 
00414   std::vector<T> _values;
00415 
00419   numeric_index_type _global_size;
00420 
00424   numeric_index_type _local_size;
00425 
00429   numeric_index_type _first_local_index;
00430 
00434   numeric_index_type _last_local_index;
00435 };
00436 
00437 
00438 //--------------------------------------------------------------------------
00439 // DistributedVector inline methods
00440 template <typename T>
00441 inline
00442 DistributedVector<T>::DistributedVector (const Parallel::Communicator &comm_in,
00443                                          const ParallelType ptype)
00444   : NumericVector<T>(comm_in, ptype),
00445     _global_size      (0),
00446     _local_size       (0),
00447     _first_local_index(0),
00448     _last_local_index (0)
00449 {
00450   this->_type = ptype;
00451 }
00452 
00453 
00454 
00455 template <typename T>
00456 inline
00457 DistributedVector<T>::DistributedVector (const Parallel::Communicator &comm_in,
00458                                          const numeric_index_type n,
00459                                          const ParallelType ptype)
00460   : NumericVector<T>(comm_in, ptype)
00461 {
00462   this->init(n, n, false, ptype);
00463 }
00464 
00465 
00466 
00467 template <typename T>
00468 inline
00469 DistributedVector<T>::DistributedVector (const Parallel::Communicator &comm_in,
00470                                          const numeric_index_type n,
00471                                          const numeric_index_type n_local,
00472                                          const ParallelType ptype)
00473   : NumericVector<T>(comm_in, ptype)
00474 {
00475   this->init(n, n_local, false, ptype);
00476 }
00477 
00478 
00479 
00480 template <typename T>
00481 inline
00482 DistributedVector<T>::DistributedVector (const Parallel::Communicator &comm_in,
00483                                          const numeric_index_type n,
00484                                          const numeric_index_type n_local,
00485                                          const std::vector<numeric_index_type>& ghost,
00486                                          const ParallelType ptype)
00487   : NumericVector<T>(comm_in, ptype)
00488 {
00489   this->init(n, n_local, ghost, false, ptype);
00490 }
00491 
00492 
00493 
00494 template <typename T>
00495 inline
00496 DistributedVector<T>::~DistributedVector ()
00497 {
00498   this->clear ();
00499 }
00500 
00501 
00502 
00503 template <typename T>
00504 inline
00505 void DistributedVector<T>::init (const numeric_index_type n,
00506                                  const numeric_index_type n_local,
00507                                  const bool fast,
00508                                  const ParallelType ptype)
00509 {
00510   // This function must be run on all processors at once
00511   parallel_object_only();
00512 
00513   libmesh_assert_less_equal (n_local, n);
00514 
00515   if (ptype == AUTOMATIC)
00516     {
00517       if (n == n_local)
00518         this->_type = SERIAL;
00519       else
00520         this->_type = PARALLEL;
00521     }
00522   else
00523     this->_type = ptype;
00524 
00525   libmesh_assert ((this->_type==SERIAL && n==n_local) ||
00526                   this->_type==PARALLEL);
00527 
00528   // Clear the data structures if already initialized
00529   if (this->initialized())
00530     this->clear();
00531 
00532   // Initialize data structures
00533   _values.resize(n_local);
00534   _local_size  = n_local;
00535   _global_size = n;
00536 
00537   _first_local_index = 0;
00538 
00539 #ifdef LIBMESH_HAVE_MPI
00540 
00541   std::vector<numeric_index_type> local_sizes (this->n_processors(), 0);
00542 
00543   local_sizes[this->processor_id()] = n_local;
00544 
00545   this->comm().sum(local_sizes);
00546 
00547   // _first_local_index is the sum of _local_size
00548   // for all processor ids less than ours
00549   for (processor_id_type p=0; p!=this->processor_id(); p++)
00550     _first_local_index += local_sizes[p];
00551 
00552 
00553 #  ifdef DEBUG
00554   // Make sure all the local sizes sum up to the global
00555   // size, otherwise there is big trouble!
00556   numeric_index_type dbg_sum=0;
00557 
00558   for (processor_id_type p=0; p!=this->n_processors(); p++)
00559     dbg_sum += local_sizes[p];
00560 
00561   libmesh_assert_equal_to (dbg_sum, n);
00562 
00563 #  endif
00564 
00565 #else
00566 
00567   // No other options without MPI!
00568   if (n != n_local)
00569     libmesh_error_msg("ERROR:  MPI is required for n != n_local!");
00570 
00571 #endif
00572 
00573   _last_local_index = _first_local_index + n_local;
00574 
00575   // Set the initialized flag
00576   this->_is_initialized = true;
00577 
00578   // Zero the components unless directed otherwise
00579   if (!fast)
00580     this->zero();
00581 }
00582 
00583 
00584 template <typename T>
00585 inline
00586 void DistributedVector<T>::init (const numeric_index_type n,
00587                                  const numeric_index_type n_local,
00588                                  const std::vector<numeric_index_type>& /*ghost*/,
00589                                  const bool fast,
00590                                  const ParallelType ptype)
00591 {
00592   // TODO: we shouldn't ignore the ghost sparsity pattern
00593   this->init(n, n_local, fast, ptype);
00594 }
00595 
00596 
00597 
00598 /* Default implementation for solver packages for which ghosted
00599    vectors are not yet implemented.  */
00600 template <class T>
00601 void DistributedVector<T>::init (const NumericVector<T>& other,
00602                                  const bool fast)
00603 {
00604   this->init(other.size(),other.local_size(),fast,other.type());
00605 }
00606 
00607 
00608 
00609 template <typename T>
00610 inline
00611 void DistributedVector<T>::init (const numeric_index_type n,
00612                                  const bool fast,
00613                                  const ParallelType ptype)
00614 {
00615   this->init(n,n,fast,ptype);
00616 }
00617 
00618 
00619 
00620 template <typename T>
00621 inline
00622 void DistributedVector<T>::close ()
00623 {
00624   libmesh_assert (this->initialized());
00625 
00626   this->_is_closed = true;
00627 }
00628 
00629 
00630 
00631 template <typename T>
00632 inline
00633 void DistributedVector<T>::clear ()
00634 {
00635   _values.clear();
00636 
00637   _global_size =
00638     _local_size =
00639     _first_local_index =
00640     _last_local_index = 0;
00641 
00642 
00643   this->_is_closed = this->_is_initialized = false;
00644 }
00645 
00646 
00647 
00648 template <typename T>
00649 inline
00650 void DistributedVector<T>::zero ()
00651 {
00652   libmesh_assert (this->initialized());
00653   libmesh_assert_equal_to (_values.size(), _local_size);
00654   libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
00655 
00656   std::fill (_values.begin(),
00657              _values.end(),
00658              0.);
00659 }
00660 
00661 
00662 
00663 template <typename T>
00664 inline
00665 UniquePtr<NumericVector<T> > DistributedVector<T>::zero_clone () const
00666 {
00667   NumericVector<T>* cloned_vector = new DistributedVector<T>(this->comm());
00668   cloned_vector->init(*this);
00669   return UniquePtr<NumericVector<T> >(cloned_vector);
00670 }
00671 
00672 
00673 
00674 template <typename T>
00675 inline
00676 UniquePtr<NumericVector<T> > DistributedVector<T>::clone () const
00677 {
00678   NumericVector<T>* cloned_vector = new DistributedVector<T>(this->comm());
00679   cloned_vector->init(*this, true);
00680   *cloned_vector = *this;
00681   return UniquePtr<NumericVector<T> >(cloned_vector);
00682 }
00683 
00684 
00685 
00686 template <typename T>
00687 inline
00688 numeric_index_type DistributedVector<T>::size () const
00689 {
00690   libmesh_assert (this->initialized());
00691   libmesh_assert_equal_to (_values.size(), _local_size);
00692   libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
00693 
00694   return _global_size;
00695 }
00696 
00697 
00698 
00699 template <typename T>
00700 inline
00701 numeric_index_type DistributedVector<T>::local_size () const
00702 {
00703   libmesh_assert (this->initialized());
00704   libmesh_assert_equal_to (_values.size(), _local_size);
00705   libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
00706 
00707   return _local_size;
00708 }
00709 
00710 
00711 
00712 template <typename T>
00713 inline
00714 numeric_index_type DistributedVector<T>::first_local_index () const
00715 {
00716   libmesh_assert (this->initialized());
00717   libmesh_assert_equal_to (_values.size(), _local_size);
00718   libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
00719 
00720   return _first_local_index;
00721 }
00722 
00723 
00724 
00725 template <typename T>
00726 inline
00727 numeric_index_type DistributedVector<T>::last_local_index () const
00728 {
00729   libmesh_assert (this->initialized());
00730   libmesh_assert_equal_to (_values.size(), _local_size);
00731   libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
00732 
00733   return _last_local_index;
00734 }
00735 
00736 
00737 
00738 template <typename T>
00739 inline
00740 T DistributedVector<T>::operator() (const numeric_index_type i) const
00741 {
00742   libmesh_assert (this->initialized());
00743   libmesh_assert_equal_to (_values.size(), _local_size);
00744   libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
00745   libmesh_assert ( ((i >= first_local_index()) &&
00746                     (i <  last_local_index())) );
00747 
00748   return _values[i - _first_local_index];
00749 }
00750 
00751 
00752 
00753 template <typename T>
00754 inline
00755 void DistributedVector<T>::set (const numeric_index_type i, const T value)
00756 {
00757   libmesh_assert (this->initialized());
00758   libmesh_assert_equal_to (_values.size(), _local_size);
00759   libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
00760   libmesh_assert_less (i, size());
00761   libmesh_assert_less (i-first_local_index(), local_size());
00762 
00763   _values[i - _first_local_index] = value;
00764 }
00765 
00766 
00767 
00768 template <typename T>
00769 inline
00770 void DistributedVector<T>::add (const numeric_index_type i, const T value)
00771 {
00772   libmesh_assert (this->initialized());
00773   libmesh_assert_equal_to (_values.size(), _local_size);
00774   libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
00775   libmesh_assert_less (i, size());
00776   libmesh_assert_less (i-first_local_index(), local_size());
00777 
00778   _values[i - _first_local_index] += value;
00779 }
00780 
00781 
00782 
00783 template <typename T>
00784 inline
00785 Real DistributedVector<T>::min () const
00786 {
00787   // This function must be run on all processors at once
00788   parallel_object_only();
00789 
00790   libmesh_assert (this->initialized());
00791   libmesh_assert_equal_to (_values.size(), _local_size);
00792   libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
00793 
00794   Real local_min = _values.size() ?
00795     libmesh_real(_values[0]) : std::numeric_limits<Real>::max();
00796   for (numeric_index_type i = 1; i < _values.size(); ++i)
00797     local_min = std::min(libmesh_real(_values[i]), local_min);
00798 
00799   this->comm().min(local_min);
00800 
00801   return local_min;
00802 }
00803 
00804 
00805 
00806 template <typename T>
00807 inline
00808 Real DistributedVector<T>::max() const
00809 {
00810   // This function must be run on all processors at once
00811   parallel_object_only();
00812 
00813   libmesh_assert (this->initialized());
00814   libmesh_assert_equal_to (_values.size(), _local_size);
00815   libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
00816 
00817   Real local_max = _values.size() ?
00818     libmesh_real(_values[0]) : -std::numeric_limits<Real>::max();
00819   for (numeric_index_type i = 1; i < _values.size(); ++i)
00820     local_max = std::max(libmesh_real(_values[i]), local_max);
00821 
00822   this->comm().max(local_max);
00823 
00824   return local_max;
00825 }
00826 
00827 
00828 template <typename T>
00829 inline
00830 void DistributedVector<T>::swap (NumericVector<T> &other)
00831 {
00832   DistributedVector<T>& v = cast_ref<DistributedVector<T>&>(other);
00833 
00834   std::swap(_global_size, v._global_size);
00835   std::swap(_local_size, v._local_size);
00836   std::swap(_first_local_index, v._first_local_index);
00837   std::swap(_last_local_index, v._last_local_index);
00838 
00839   // This should be O(1) with any reasonable STL implementation
00840   std::swap(_values, v._values);
00841 }
00842 
00843 } // namespace libMesh
00844 
00845 
00846 #endif  // LIBMESH_DISTRIBUTED_VECTOR_H