$extrastylesheet
laspack_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 
00021 #ifndef LIBMESH_LASPACK_VECTOR_H
00022 #define LIBMESH_LASPACK_VECTOR_H
00023 
00024 
00025 
00026 #include "libmesh/libmesh_common.h"
00027 
00028 #ifdef LIBMESH_HAVE_LASPACK
00029 
00030 // Local includes
00031 #include "libmesh/numeric_vector.h"
00032 
00033 // C++ includes
00034 #include <cstdio> // for std::sprintf
00035 
00036 // Laspack includes
00037 #include <operats.h>
00038 #include <qvector.h>
00039 
00040 namespace libMesh
00041 {
00042 
00043 
00044 // Forward declarations
00045 template <typename T> class LaspackLinearSolver;
00046 template <typename T> class SparseMatrix;
00047 
00048 
00056 template <typename T>
00057 class LaspackVector : public NumericVector<T>
00058 {
00059 public:
00060 
00064   explicit
00065   LaspackVector (const Parallel::Communicator &comm,
00066                  const ParallelType = AUTOMATIC);
00067 
00071   explicit
00072   LaspackVector (const Parallel::Communicator &comm,
00073                  const numeric_index_type n,
00074                  const ParallelType = AUTOMATIC);
00075 
00080   LaspackVector (const Parallel::Communicator &comm,
00081                  const numeric_index_type n,
00082                  const numeric_index_type n_local,
00083                  const ParallelType = AUTOMATIC);
00084 
00090   LaspackVector (const Parallel::Communicator &comm,
00091                  const numeric_index_type N,
00092                  const numeric_index_type n_local,
00093                  const std::vector<numeric_index_type>& ghost,
00094                  const ParallelType = AUTOMATIC);
00095 
00100   ~LaspackVector ();
00101 
00105   void close ();
00106 
00110   void clear ();
00111 
00116   void zero ();
00117 
00123   virtual UniquePtr<NumericVector<T> > zero_clone () const;
00124 
00128   UniquePtr<NumericVector<T> > clone () const;
00129 
00143   void init (const numeric_index_type N,
00144              const numeric_index_type n_local,
00145              const bool         fast=false,
00146              const ParallelType ptype=AUTOMATIC);
00147 
00151   void init (const numeric_index_type N,
00152              const bool         fast=false,
00153              const ParallelType ptype=AUTOMATIC);
00154 
00159   void init (const numeric_index_type /*N*/,
00160              const numeric_index_type /*n_local*/,
00161              const std::vector<numeric_index_type>& /*ghost*/,
00162              const bool /*fast*/ = false,
00163              const ParallelType = AUTOMATIC);
00164 
00169   virtual void init (const NumericVector<T>& other,
00170                      const bool fast = false);
00171 
00175   NumericVector<T> & operator= (const T s);
00176 
00180   NumericVector<T> & operator= (const NumericVector<T> &V);
00181 
00185   LaspackVector<T> & operator= (const LaspackVector<T> &V);
00186 
00190   NumericVector<T> & operator= (const std::vector<T> &v);
00191 
00197   Real min () const;
00198 
00204   Real max () const;
00205 
00209   T sum () const;
00210 
00215   Real l1_norm () const;
00216 
00222   Real l2_norm () const;
00223 
00229   Real linfty_norm () const;
00230 
00238   numeric_index_type size () const;
00239 
00244   numeric_index_type local_size() const;
00245 
00250   numeric_index_type first_local_index() const;
00251 
00256   numeric_index_type last_local_index() const;
00257 
00261   T operator() (const numeric_index_type i) const;
00262 
00267   NumericVector<T> & operator += (const NumericVector<T> &V);
00268 
00273   NumericVector<T> & operator -= (const NumericVector<T> &V);
00274 
00278   virtual NumericVector<T> & operator /= (NumericVector<T> & v);
00279 
00283   virtual void reciprocal();
00284 
00289   virtual void conjugate();
00290 
00294   void set (const numeric_index_type i, const T value);
00295 
00299   void add (const numeric_index_type i, const T value);
00300 
00306   void add (const T s);
00307 
00313   void add (const NumericVector<T>& V);
00314 
00320   void add (const T a, const NumericVector<T>& v);
00321 
00326   using NumericVector<T>::add_vector;
00327 
00332   void add_vector (const NumericVector<T> &,
00333                    const SparseMatrix<T> &);
00334 
00339   void add_vector_transpose (const NumericVector<T> &,
00340                              const SparseMatrix<T> &);
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 
00387 
00394   void localize_to_one (std::vector<T>& v_local,
00395                         const processor_id_type proc_id=0) const;
00396 
00401   virtual void pointwise_mult (const NumericVector<T>& vec1,
00402                                const NumericVector<T>& vec2);
00403 
00407   virtual void swap (NumericVector<T> &v);
00408 
00409 private:
00410 
00415   QVector _vec;
00416 
00420   friend class LaspackLinearSolver<T>;
00421 };
00422 
00423 
00424 
00425 //----------------------- ----------------------------------
00426 // LaspackVector inline methods
00427 template <typename T>
00428 inline
00429 LaspackVector<T>::LaspackVector (const Parallel::Communicator &comm,
00430                                  const ParallelType ptype)
00431   : NumericVector<T>(comm, ptype)
00432 {
00433   this->_type = ptype;
00434 }
00435 
00436 
00437 
00438 template <typename T>
00439 inline
00440 LaspackVector<T>::LaspackVector (const Parallel::Communicator &comm,
00441                                  const numeric_index_type n,
00442                                  const ParallelType ptype)
00443   : NumericVector<T>(comm, ptype)
00444 {
00445   this->init(n, n, false, ptype);
00446 }
00447 
00448 
00449 
00450 template <typename T>
00451 inline
00452 LaspackVector<T>::LaspackVector (const Parallel::Communicator &comm,
00453                                  const numeric_index_type n,
00454                                  const numeric_index_type n_local,
00455                                  const ParallelType ptype)
00456   : NumericVector<T>(comm, ptype)
00457 {
00458   this->init(n, n_local, false, ptype);
00459 }
00460 
00461 
00462 
00463 template <typename T>
00464 inline
00465 LaspackVector<T>::LaspackVector (const Parallel::Communicator &comm,
00466                                  const numeric_index_type N,
00467                                  const numeric_index_type n_local,
00468                                  const std::vector<numeric_index_type>& ghost,
00469                                  const ParallelType ptype)
00470   : NumericVector<T>(comm, ptype)
00471 {
00472   this->init(N, n_local, ghost, false, ptype);
00473 }
00474 
00475 
00476 
00477 template <typename T>
00478 inline
00479 LaspackVector<T>::~LaspackVector ()
00480 {
00481   this->clear ();
00482 }
00483 
00484 
00485 
00486 template <typename T>
00487 inline
00488 void LaspackVector<T>::init (const numeric_index_type n,
00489                              const numeric_index_type libmesh_dbg_var(n_local),
00490                              const bool fast,
00491                              const ParallelType)
00492 {
00493   // Laspack vectors only for serial cases,
00494   // but can provide a "parallel" vector on one processor.
00495   libmesh_assert_equal_to (n, n_local);
00496 
00497   this->_type = SERIAL;
00498 
00499   // Clear initialized vectors
00500   if (this->initialized())
00501     this->clear();
00502 
00503   // create a sequential vector
00504 
00505   static int cnt = 0;
00506   char foo[80];
00507   std::sprintf(foo,  "Vec-%d", cnt++);
00508 
00509   V_Constr(&_vec, const_cast<char*>(foo), n, Normal, _LPTrue);
00510 
00511   this->_is_initialized = true;
00512 #ifndef NDEBUG
00513   this->_is_closed = true;
00514 #endif
00515 
00516   // Optionally zero out all components
00517   if (fast == false)
00518     this->zero ();
00519 
00520   return;
00521 }
00522 
00523 
00524 
00525 template <typename T>
00526 inline
00527 void LaspackVector<T>::init (const numeric_index_type n,
00528                              const bool fast,
00529                              const ParallelType ptype)
00530 {
00531   this->init(n,n,fast,ptype);
00532 }
00533 
00534 
00535 template <typename T>
00536 inline
00537 void LaspackVector<T>::init (const numeric_index_type n,
00538                              const numeric_index_type n_local,
00539                              const std::vector<numeric_index_type>& libmesh_dbg_var(ghost),
00540                              const bool fast,
00541                              const ParallelType ptype)
00542 {
00543   libmesh_assert(ghost.empty());
00544   this->init(n,n_local,fast,ptype);
00545 }
00546 
00547 
00548 
00549 /* Default implementation for solver packages for which ghosted
00550    vectors are not yet implemented.  */
00551 template <class T>
00552 void LaspackVector<T>::init (const NumericVector<T>& other,
00553                              const bool fast)
00554 {
00555   this->init(other.size(),other.local_size(),fast,other.type());
00556 }
00557 
00558 
00559 
00560 template <typename T>
00561 inline
00562 void LaspackVector<T>::close ()
00563 {
00564   libmesh_assert (this->initialized());
00565 
00566 #ifndef NDEBUG
00567   this->_is_closed = true;
00568 #endif
00569 }
00570 
00571 
00572 
00573 template <typename T>
00574 inline
00575 void LaspackVector<T>::clear ()
00576 {
00577   if (this->initialized())
00578     {
00579       V_Destr (&_vec);
00580     }
00581 
00582   this->_is_initialized = false;
00583 #ifndef NDEBUG
00584   this->_is_closed = false;
00585 #endif
00586 }
00587 
00588 
00589 
00590 template <typename T> inline
00591 void LaspackVector<T>::zero ()
00592 {
00593   libmesh_assert (this->initialized());
00594   libmesh_assert (this->closed());
00595 
00596   V_SetAllCmp (&_vec, 0.);
00597 }
00598 
00599 
00600 
00601 template <typename T>
00602 inline
00603 UniquePtr<NumericVector<T> > LaspackVector<T>::zero_clone () const
00604 {
00605   NumericVector<T> * cloned_vector = new LaspackVector<T>(this->comm());
00606 
00607   cloned_vector->init(*this);
00608 
00609   return UniquePtr<NumericVector<T> >(cloned_vector);
00610 }
00611 
00612 
00613 
00614 template <typename T>
00615 inline
00616 UniquePtr<NumericVector<T> > LaspackVector<T>::clone () const
00617 {
00618   NumericVector<T> * cloned_vector = new LaspackVector<T>(this->comm());
00619 
00620   cloned_vector->init(*this, true);
00621 
00622   *cloned_vector = *this;
00623 
00624   return UniquePtr<NumericVector<T> >(cloned_vector);
00625 }
00626 
00627 
00628 
00629 template <typename T>
00630 inline
00631 numeric_index_type LaspackVector<T>::size () const
00632 {
00633   libmesh_assert (this->initialized());
00634 
00635   return static_cast<numeric_index_type>(V_GetDim(const_cast<QVector*>(&_vec)));
00636 }
00637 
00638 
00639 
00640 template <typename T>
00641 inline
00642 numeric_index_type LaspackVector<T>::local_size () const
00643 {
00644   libmesh_assert (this->initialized());
00645 
00646   return this->size();
00647 }
00648 
00649 
00650 
00651 template <typename T>
00652 inline
00653 numeric_index_type LaspackVector<T>::first_local_index () const
00654 {
00655   libmesh_assert (this->initialized());
00656 
00657   return 0;
00658 }
00659 
00660 
00661 
00662 template <typename T>
00663 inline
00664 numeric_index_type LaspackVector<T>::last_local_index () const
00665 {
00666   libmesh_assert (this->initialized());
00667 
00668   return this->size();
00669 }
00670 
00671 
00672 
00673 template <typename T>
00674 inline
00675 void LaspackVector<T>::set (const numeric_index_type i, const T value)
00676 {
00677   libmesh_assert (this->initialized());
00678   libmesh_assert_less (i, this->size());
00679 
00680   V_SetCmp (&_vec, i+1, value);
00681 
00682 #ifndef NDEBUG
00683   this->_is_closed = false;
00684 #endif
00685 }
00686 
00687 
00688 
00689 template <typename T>
00690 inline
00691 void LaspackVector<T>::add (const numeric_index_type i, const T value)
00692 {
00693   libmesh_assert (this->initialized());
00694   libmesh_assert_less (i, this->size());
00695 
00696   V_AddCmp (&_vec, i+1, value);
00697 
00698 #ifndef NDEBUG
00699   this->_is_closed = false;
00700 #endif
00701 }
00702 
00703 
00704 
00705 template <typename T>
00706 inline
00707 T LaspackVector<T>::operator() (const numeric_index_type i) const
00708 {
00709   libmesh_assert (this->initialized());
00710   libmesh_assert ( ((i >= this->first_local_index()) &&
00711                     (i <  this->last_local_index())) );
00712 
00713 
00714   return static_cast<T>(V_GetCmp(const_cast<QVector*>(&_vec), i+1));
00715 }
00716 
00717 
00718 
00719 template <typename T>
00720 inline
00721 void LaspackVector<T>::swap (NumericVector<T> &other)
00722 {
00723   LaspackVector<T>& v = cast_ref<LaspackVector<T>&>(other);
00724 
00725   // This is all grossly dependent on Laspack version...
00726 
00727   std::swap(_vec.Name, v._vec.Name);
00728   std::swap(_vec.Dim, v._vec.Dim);
00729   std::swap(_vec.Instance, v._vec.Instance);
00730   std::swap(_vec.LockLevel, v._vec.LockLevel);
00731   std::swap(_vec.Multipl, v._vec.Multipl);
00732   std::swap(_vec.OwnData, v._vec.OwnData);
00733 
00734   // This should still be O(1), since _vec.Cmp is just a pointer to
00735   // data on the heap
00736 
00737   std::swap(_vec.Cmp, v._vec.Cmp);
00738 }
00739 
00740 
00741 } // namespace libMesh
00742 
00743 
00744 #endif // #ifdef LIBMESH_HAVE_LASPACK
00745 #endif // LIBMESH_LASPACK_VECTOR_H