$extrastylesheet
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_PETSC_VECTOR_H 00022 #define LIBMESH_PETSC_VECTOR_H 00023 00024 00025 #include "libmesh/libmesh_config.h" 00026 00027 #ifdef LIBMESH_HAVE_PETSC 00028 00029 // Local includes 00030 #include "libmesh/numeric_vector.h" 00031 #include "libmesh/petsc_macro.h" 00032 #include "libmesh/libmesh_common.h" 00033 #include LIBMESH_INCLUDE_UNORDERED_MAP 00034 00038 EXTERN_C_FOR_PETSC_BEGIN 00039 # include <petscvec.h> 00040 EXTERN_C_FOR_PETSC_END 00041 00042 // C++ includes 00043 #include <cstddef> 00044 #include <cstring> 00045 #include <vector> 00046 00047 namespace libMesh 00048 { 00049 00050 00051 00052 // forward declarations 00053 template <typename T> class SparseMatrix; 00054 00062 template <typename T> 00063 class PetscVector : public NumericVector<T> 00064 { 00065 public: 00066 00070 explicit 00071 PetscVector (const Parallel::Communicator &comm_in, 00072 const ParallelType type = AUTOMATIC); 00073 00077 explicit 00078 PetscVector (const Parallel::Communicator &comm_in, 00079 const numeric_index_type n, 00080 const ParallelType type = AUTOMATIC); 00081 00086 PetscVector (const Parallel::Communicator &comm_in, 00087 const numeric_index_type n, 00088 const numeric_index_type n_local, 00089 const ParallelType type = AUTOMATIC); 00090 00096 PetscVector (const Parallel::Communicator &comm_in, 00097 const numeric_index_type N, 00098 const numeric_index_type n_local, 00099 const std::vector<numeric_index_type>& ghost, 00100 const ParallelType type = AUTOMATIC); 00101 00109 PetscVector(Vec v, 00110 const Parallel::Communicator &comm_in 00111 LIBMESH_CAN_DEFAULT_TO_COMMWORLD); 00112 00117 ~PetscVector (); 00118 00122 void close (); 00123 00127 void clear (); 00128 00133 void zero (); 00134 00140 virtual UniquePtr<NumericVector<T> > zero_clone () const; 00141 00145 UniquePtr<NumericVector<T> > clone () const; 00146 00160 void init (const numeric_index_type N, 00161 const numeric_index_type n_local, 00162 const bool fast=false, 00163 const ParallelType type=AUTOMATIC); 00164 00168 void init (const numeric_index_type N, 00169 const bool fast=false, 00170 const ParallelType type=AUTOMATIC); 00171 00176 virtual void init (const numeric_index_type /*N*/, 00177 const numeric_index_type /*n_local*/, 00178 const std::vector<numeric_index_type>& /*ghost*/, 00179 const bool /*fast*/ = false, 00180 const ParallelType = AUTOMATIC); 00181 00186 virtual void init (const NumericVector<T>& other, 00187 const bool fast = false); 00188 00189 // /** 00190 // * Change the dimension to that of the 00191 // * vector \p V. The same applies as for 00192 // * the other \p init function. 00193 // * 00194 // * The elements of \p V are not copied, i.e. 00195 // * this function is the same as calling 00196 // * \p init(V.size(),fast). 00197 // */ 00198 // void init (const NumericVector<T>& V, 00199 // const bool fast=false); 00200 00204 NumericVector<T> & operator= (const T s); 00205 00209 NumericVector<T> & operator= (const NumericVector<T> &V); 00210 00214 PetscVector<T> & operator= (const PetscVector<T> &V); 00215 00219 NumericVector<T> & operator= (const std::vector<T> &v); 00220 00226 Real min () const; 00227 00233 Real max () const; 00234 00238 T sum () const; 00239 00244 Real l1_norm () const; 00245 00251 Real l2_norm () const; 00252 00258 Real linfty_norm () const; 00259 00267 numeric_index_type size () const; 00268 00273 numeric_index_type local_size() const; 00274 00279 numeric_index_type first_local_index() const; 00280 00285 numeric_index_type last_local_index() const; 00286 00293 numeric_index_type map_global_to_local_index(const numeric_index_type i) const; 00294 00298 T operator() (const numeric_index_type i) const; 00299 00306 virtual void get(const std::vector<numeric_index_type>& index, T* values) const; 00307 00312 NumericVector<T> & operator += (const NumericVector<T> &V); 00313 00318 NumericVector<T> & operator -= (const NumericVector<T> &V); 00319 00323 virtual void reciprocal(); 00324 00329 virtual void conjugate(); 00330 00334 void set (const numeric_index_type i, const T value); 00335 00339 void add (const numeric_index_type i, const T value); 00340 00346 void add (const T s); 00347 00353 void add (const NumericVector<T>& V); 00354 00360 void add (const T a, const NumericVector<T>& v); 00361 00366 using NumericVector<T>::add_vector; 00367 00372 void add_vector (const T* v, 00373 const std::vector<numeric_index_type>& dof_indices); 00374 00379 void add_vector (const NumericVector<T> &V, 00380 const SparseMatrix<T> &A); 00381 00387 void add_vector_transpose (const NumericVector<T> &V, 00388 const SparseMatrix<T> &A_trans); 00389 00395 void add_vector_conjugate_transpose (const NumericVector<T> &V, 00396 const SparseMatrix<T> &A_trans); 00397 00402 using NumericVector<T>::insert; 00403 00408 virtual void insert (const T* v, 00409 const std::vector<numeric_index_type>& dof_indices); 00410 00415 void scale (const T factor); 00416 00420 virtual NumericVector<T> & operator /= (NumericVector<T> & v); 00421 00426 virtual void abs(); 00427 00432 virtual T dot(const NumericVector<T>&) const; 00433 00438 virtual T indefinite_dot(const NumericVector<T>&) const; 00439 00444 void localize (std::vector<T>& v_local) const; 00445 00450 void localize (NumericVector<T>& v_local) const; 00451 00457 void localize (NumericVector<T>& v_local, 00458 const std::vector<numeric_index_type>& send_list) const; 00459 00464 void localize (const numeric_index_type first_local_idx, 00465 const numeric_index_type last_local_idx, 00466 const std::vector<numeric_index_type>& send_list); 00467 00474 void localize_to_one (std::vector<T>& v_local, 00475 const processor_id_type proc_id=0) const; 00476 00481 virtual void pointwise_mult (const NumericVector<T>& vec1, 00482 const NumericVector<T>& vec2); 00483 00490 virtual void print_matlab(const std::string& name = "") const; 00491 00496 virtual void create_subvector(NumericVector<T>& subvector, 00497 const std::vector<numeric_index_type>& rows) const; 00498 00502 virtual void swap (NumericVector<T> &v); 00503 00509 Vec vec () { libmesh_assert (_vec); return _vec; } 00510 00511 00512 00513 private: 00514 00519 Vec _vec; 00520 00526 mutable bool _array_is_present; 00527 00533 mutable numeric_index_type _first; 00534 00540 mutable numeric_index_type _last; 00541 00542 #ifndef NDEBUG 00543 00548 mutable numeric_index_type _local_size; 00549 #endif 00550 00556 mutable Vec _local_form; 00557 00562 mutable PetscScalar* _values; 00563 00568 void _get_array(void) const; 00569 00574 void _restore_array(void) const; 00575 00579 typedef LIBMESH_BEST_UNORDERED_MAP<numeric_index_type,numeric_index_type> GlobalToLocalMap; 00580 00585 GlobalToLocalMap _global_to_local_map; 00586 00591 bool _destroy_vec_on_exit; 00592 }; 00593 00594 00595 /*----------------------- Inline functions ----------------------------------*/ 00596 00597 00598 00599 template <typename T> 00600 inline 00601 PetscVector<T>::PetscVector (const Parallel::Communicator &comm_in, const ParallelType ptype) 00602 : NumericVector<T>(comm_in, ptype), 00603 _array_is_present(false), 00604 _first(0), 00605 _last(0), 00606 _local_form(NULL), 00607 _values(NULL), 00608 _global_to_local_map(), 00609 _destroy_vec_on_exit(true) 00610 { 00611 this->_type = ptype; 00612 } 00613 00614 00615 00616 template <typename T> 00617 inline 00618 PetscVector<T>::PetscVector (const Parallel::Communicator &comm_in, 00619 const numeric_index_type n, 00620 const ParallelType ptype) 00621 : NumericVector<T>(comm_in, ptype), 00622 _array_is_present(false), 00623 _local_form(NULL), 00624 _values(NULL), 00625 _global_to_local_map(), 00626 _destroy_vec_on_exit(true) 00627 { 00628 this->init(n, n, false, ptype); 00629 } 00630 00631 00632 00633 template <typename T> 00634 inline 00635 PetscVector<T>::PetscVector (const Parallel::Communicator &comm_in, 00636 const numeric_index_type n, 00637 const numeric_index_type n_local, 00638 const ParallelType ptype) 00639 : NumericVector<T>(comm_in, ptype), 00640 _array_is_present(false), 00641 _local_form(NULL), 00642 _values(NULL), 00643 _global_to_local_map(), 00644 _destroy_vec_on_exit(true) 00645 { 00646 this->init(n, n_local, false, ptype); 00647 } 00648 00649 00650 00651 template <typename T> 00652 inline 00653 PetscVector<T>::PetscVector (const Parallel::Communicator &comm_in, 00654 const numeric_index_type n, 00655 const numeric_index_type n_local, 00656 const std::vector<numeric_index_type>& ghost, 00657 const ParallelType ptype) 00658 : NumericVector<T>(comm_in, ptype), 00659 _array_is_present(false), 00660 _local_form(NULL), 00661 _values(NULL), 00662 _global_to_local_map(), 00663 _destroy_vec_on_exit(true) 00664 { 00665 this->init(n, n_local, ghost, false, ptype); 00666 } 00667 00668 00669 00670 00671 00672 template <typename T> 00673 inline 00674 PetscVector<T>::PetscVector (Vec v, 00675 const Parallel::Communicator &comm_in) 00676 : NumericVector<T>(comm_in, AUTOMATIC), 00677 _array_is_present(false), 00678 _local_form(NULL), 00679 _values(NULL), 00680 _global_to_local_map(), 00681 _destroy_vec_on_exit(false) 00682 { 00683 this->_vec = v; 00684 this->_is_closed = true; 00685 this->_is_initialized = true; 00686 00687 /* We need to ask PETSc about the (local to global) ghost value 00688 mapping and create the inverse mapping out of it. */ 00689 PetscErrorCode ierr=0; 00690 PetscInt petsc_local_size=0; 00691 ierr = VecGetLocalSize(_vec, &petsc_local_size); 00692 LIBMESH_CHKERRABORT(ierr); 00693 00694 // Get the vector type from PETSc. 00695 // As of Petsc 3.0.0, the VecType #define lost its const-ness, so we 00696 // need to have it in the code 00697 #if PETSC_VERSION_LESS_THAN(3,0,0) || !PETSC_VERSION_LESS_THAN(3,4,0) 00698 // Pre-3.0 and petsc-dev (as of October 2012) use non-const versions 00699 VecType ptype; 00700 #else 00701 const VecType ptype; 00702 #endif 00703 ierr = VecGetType(_vec, &ptype); 00704 LIBMESH_CHKERRABORT(ierr); 00705 00706 if((std::strcmp(ptype,VECSHARED) == 0) || (std::strcmp(ptype,VECMPI) == 0)) 00707 { 00708 #if PETSC_RELEASE_LESS_THAN(3,1,1) 00709 ISLocalToGlobalMapping mapping = _vec->mapping; 00710 #else 00711 ISLocalToGlobalMapping mapping; 00712 ierr = VecGetLocalToGlobalMapping(_vec, &mapping); 00713 LIBMESH_CHKERRABORT(ierr); 00714 #endif 00715 00716 // If is a sparsely stored vector, set up our new mapping 00717 if (mapping) 00718 { 00719 const numeric_index_type my_local_size = static_cast<numeric_index_type>(petsc_local_size); 00720 const numeric_index_type ghost_begin = static_cast<numeric_index_type>(petsc_local_size); 00721 #if PETSC_RELEASE_LESS_THAN(3,4,0) 00722 const numeric_index_type ghost_end = static_cast<numeric_index_type>(mapping->n); 00723 #else 00724 PetscInt n; 00725 ierr = ISLocalToGlobalMappingGetSize(mapping, &n); 00726 LIBMESH_CHKERRABORT(ierr); 00727 const numeric_index_type ghost_end = static_cast<numeric_index_type>(n); 00728 #endif 00729 #if PETSC_RELEASE_LESS_THAN(3,1,1) 00730 const PetscInt *indices = mapping->indices; 00731 #else 00732 const PetscInt *indices; 00733 ierr = ISLocalToGlobalMappingGetIndices(mapping,&indices); 00734 LIBMESH_CHKERRABORT(ierr); 00735 #endif 00736 for(numeric_index_type i=ghost_begin; i<ghost_end; i++) 00737 _global_to_local_map[indices[i]] = i-my_local_size; 00738 this->_type = GHOSTED; 00739 #if !PETSC_RELEASE_LESS_THAN(3,1,1) 00740 ierr = ISLocalToGlobalMappingRestoreIndices(mapping, &indices); 00741 LIBMESH_CHKERRABORT(ierr); 00742 #endif 00743 } 00744 else 00745 this->_type = PARALLEL; 00746 } 00747 else 00748 this->_type = SERIAL; 00749 00750 this->close(); 00751 } 00752 00753 00754 00755 00756 template <typename T> 00757 inline 00758 PetscVector<T>::~PetscVector () 00759 { 00760 this->clear (); 00761 } 00762 00763 00764 00765 template <typename T> 00766 inline 00767 void PetscVector<T>::init (const numeric_index_type n, 00768 const numeric_index_type n_local, 00769 const bool fast, 00770 const ParallelType ptype) 00771 { 00772 PetscErrorCode ierr=0; 00773 PetscInt petsc_n=static_cast<PetscInt>(n); 00774 PetscInt petsc_n_local=static_cast<PetscInt>(n_local); 00775 00776 00777 // Clear initialized vectors 00778 if (this->initialized()) 00779 this->clear(); 00780 00781 if (ptype == AUTOMATIC) 00782 { 00783 if (n == n_local) 00784 this->_type = SERIAL; 00785 else 00786 this->_type = PARALLEL; 00787 } 00788 else 00789 this->_type = ptype; 00790 00791 libmesh_assert ((this->_type==SERIAL && n==n_local) || 00792 this->_type==PARALLEL); 00793 00794 // create a sequential vector if on only 1 processor 00795 if (this->_type == SERIAL) 00796 { 00797 ierr = VecCreateSeq (PETSC_COMM_SELF, petsc_n, &_vec); 00798 CHKERRABORT(PETSC_COMM_SELF,ierr); 00799 00800 ierr = VecSetFromOptions (_vec); 00801 CHKERRABORT(PETSC_COMM_SELF,ierr); 00802 } 00803 // otherwise create an MPI-enabled vector 00804 else if (this->_type == PARALLEL) 00805 { 00806 #ifdef LIBMESH_HAVE_MPI 00807 libmesh_assert_less_equal (n_local, n); 00808 ierr = VecCreateMPI (this->comm().get(), petsc_n_local, petsc_n, 00809 &_vec); 00810 LIBMESH_CHKERRABORT(ierr); 00811 #else 00812 libmesh_assert_equal_to (n_local, n); 00813 ierr = VecCreateSeq (PETSC_COMM_SELF, petsc_n, &_vec); 00814 CHKERRABORT(PETSC_COMM_SELF,ierr); 00815 #endif 00816 00817 ierr = VecSetFromOptions (_vec); 00818 LIBMESH_CHKERRABORT(ierr); 00819 } 00820 else 00821 libmesh_error_msg("Unsupported type " << this->_type); 00822 00823 this->_is_initialized = true; 00824 this->_is_closed = true; 00825 00826 00827 if (fast == false) 00828 this->zero (); 00829 } 00830 00831 00832 00833 template <typename T> 00834 inline 00835 void PetscVector<T>::init (const numeric_index_type n, 00836 const bool fast, 00837 const ParallelType ptype) 00838 { 00839 this->init(n,n,fast,ptype); 00840 } 00841 00842 00843 00844 template <typename T> 00845 inline 00846 void PetscVector<T>::init (const numeric_index_type n, 00847 const numeric_index_type n_local, 00848 const std::vector<numeric_index_type>& ghost, 00849 const bool fast, 00850 const ParallelType libmesh_dbg_var(ptype)) 00851 { 00852 PetscErrorCode ierr=0; 00853 PetscInt petsc_n=static_cast<PetscInt>(n); 00854 PetscInt petsc_n_local=static_cast<PetscInt>(n_local); 00855 PetscInt petsc_n_ghost=static_cast<PetscInt>(ghost.size()); 00856 00857 // If the mesh is not disjoint, every processor will either have 00858 // all the dofs, none of the dofs, or some non-zero dofs at the 00859 // boundary between processors. 00860 // 00861 // However we can't assert this, because someone might want to 00862 // construct a GHOSTED vector which doesn't include neighbor element 00863 // dofs. Boyce tried to do so in user code, and we're going to want 00864 // to do so in System::project_vector(). 00865 // 00866 // libmesh_assert(n_local == 0 || n_local == n || !ghost.empty()); 00867 00868 libmesh_assert(sizeof(PetscInt) == sizeof(numeric_index_type)); 00869 00870 PetscInt* petsc_ghost = ghost.empty() ? PETSC_NULL : 00871 const_cast<PetscInt*>(reinterpret_cast<const PetscInt*>(&ghost[0])); 00872 00873 // Clear initialized vectors 00874 if (this->initialized()) 00875 this->clear(); 00876 00877 libmesh_assert(ptype == AUTOMATIC || ptype == GHOSTED); 00878 this->_type = GHOSTED; 00879 00880 /* Make the global-to-local ghost cell map. */ 00881 for(numeric_index_type i=0; i<ghost.size(); i++) 00882 { 00883 _global_to_local_map[ghost[i]] = i; 00884 } 00885 00886 /* Create vector. */ 00887 ierr = VecCreateGhost (this->comm().get(), petsc_n_local, petsc_n, 00888 petsc_n_ghost, petsc_ghost, 00889 &_vec); 00890 LIBMESH_CHKERRABORT(ierr); 00891 00892 ierr = VecSetFromOptions (_vec); 00893 LIBMESH_CHKERRABORT(ierr); 00894 00895 this->_is_initialized = true; 00896 this->_is_closed = true; 00897 00898 if (fast == false) 00899 this->zero (); 00900 } 00901 00902 00903 00904 template <typename T> 00905 inline 00906 void PetscVector<T>::init (const NumericVector<T>& other, 00907 const bool fast) 00908 { 00909 // Clear initialized vectors 00910 if (this->initialized()) 00911 this->clear(); 00912 00913 const PetscVector<T>& v = cast_ref<const PetscVector<T>&>(other); 00914 00915 // Other vector should restore array. 00916 if(v.initialized()) 00917 { 00918 v._restore_array(); 00919 } 00920 00921 this->_global_to_local_map = v._global_to_local_map; 00922 00923 // Even if we're initializeing sizes based on an uninitialized or 00924 // unclosed vector, *this* vector is being initialized now and is 00925 // initially closed. 00926 this->_is_closed = true; // v._is_closed; 00927 this->_is_initialized = true; // v._is_initialized; 00928 00929 this->_type = v._type; 00930 00931 // We want to have a valid Vec, even if it's initially of size zero 00932 00933 // if (v.size() != 0) 00934 // { 00935 PetscErrorCode ierr = 0; 00936 00937 ierr = VecDuplicate (v._vec, &this->_vec); 00938 LIBMESH_CHKERRABORT(ierr); 00939 // } 00940 00941 if (fast == false) 00942 this->zero (); 00943 } 00944 00945 00946 00947 template <typename T> 00948 inline 00949 void PetscVector<T>::close () 00950 { 00951 this->_restore_array(); 00952 00953 PetscErrorCode ierr=0; 00954 00955 ierr = VecAssemblyBegin(_vec); 00956 LIBMESH_CHKERRABORT(ierr); 00957 ierr = VecAssemblyEnd(_vec); 00958 LIBMESH_CHKERRABORT(ierr); 00959 00960 if(this->type() == GHOSTED) 00961 { 00962 ierr = VecGhostUpdateBegin(_vec,INSERT_VALUES,SCATTER_FORWARD); 00963 LIBMESH_CHKERRABORT(ierr); 00964 ierr = VecGhostUpdateEnd(_vec,INSERT_VALUES,SCATTER_FORWARD); 00965 LIBMESH_CHKERRABORT(ierr); 00966 } 00967 00968 this->_is_closed = true; 00969 } 00970 00971 00972 00973 template <typename T> 00974 inline 00975 void PetscVector<T>::clear () 00976 { 00977 if (this->initialized()) 00978 this->_restore_array(); 00979 00980 if ((this->initialized()) && (this->_destroy_vec_on_exit)) 00981 { 00982 PetscErrorCode ierr=0; 00983 00984 ierr = LibMeshVecDestroy(&_vec); 00985 LIBMESH_CHKERRABORT(ierr); 00986 } 00987 00988 this->_is_closed = this->_is_initialized = false; 00989 00990 _global_to_local_map.clear(); 00991 } 00992 00993 00994 00995 template <typename T> 00996 inline 00997 void PetscVector<T>::zero () 00998 { 00999 libmesh_assert(this->closed()); 01000 01001 this->_restore_array(); 01002 01003 PetscErrorCode ierr=0; 01004 01005 PetscScalar z=0.; 01006 01007 if(this->type() != GHOSTED) 01008 { 01009 #if PETSC_VERSION_LESS_THAN(2,3,0) 01010 // 2.2.x & earlier style 01011 ierr = VecSet (&z, _vec); 01012 LIBMESH_CHKERRABORT(ierr); 01013 #else 01014 // 2.3.x & newer 01015 ierr = VecSet (_vec, z); 01016 LIBMESH_CHKERRABORT(ierr); 01017 #endif 01018 } 01019 else 01020 { 01021 /* Vectors that include ghost values require a special 01022 handling. */ 01023 Vec loc_vec; 01024 ierr = VecGhostGetLocalForm (_vec,&loc_vec); 01025 LIBMESH_CHKERRABORT(ierr); 01026 #if PETSC_VERSION_LESS_THAN(2,3,0) 01027 // 2.2.x & earlier style 01028 ierr = VecSet (&z, loc_vec); 01029 LIBMESH_CHKERRABORT(ierr); 01030 #else 01031 // 2.3.x & newer 01032 ierr = VecSet (loc_vec, z); 01033 LIBMESH_CHKERRABORT(ierr); 01034 #endif 01035 ierr = VecGhostRestoreLocalForm (_vec,&loc_vec); 01036 LIBMESH_CHKERRABORT(ierr); 01037 } 01038 } 01039 01040 01041 01042 template <typename T> 01043 inline 01044 UniquePtr<NumericVector<T> > PetscVector<T>::zero_clone () const 01045 { 01046 NumericVector<T>* cloned_vector = new PetscVector<T>(this->comm(), this->type()); 01047 cloned_vector->init(*this); 01048 return UniquePtr<NumericVector<T> >(cloned_vector); 01049 } 01050 01051 01052 01053 template <typename T> 01054 inline 01055 UniquePtr<NumericVector<T> > PetscVector<T>::clone () const 01056 { 01057 NumericVector<T>* cloned_vector = new PetscVector<T>(this->comm(), this->type()); 01058 cloned_vector->init(*this, true); 01059 *cloned_vector = *this; 01060 return UniquePtr<NumericVector<T> >(cloned_vector); 01061 } 01062 01063 01064 01065 template <typename T> 01066 inline 01067 numeric_index_type PetscVector<T>::size () const 01068 { 01069 libmesh_assert (this->initialized()); 01070 01071 PetscErrorCode ierr=0; 01072 PetscInt petsc_size=0; 01073 01074 if (!this->initialized()) 01075 return 0; 01076 01077 ierr = VecGetSize(_vec, &petsc_size); 01078 LIBMESH_CHKERRABORT(ierr); 01079 01080 return static_cast<numeric_index_type>(petsc_size); 01081 } 01082 01083 01084 01085 template <typename T> 01086 inline 01087 numeric_index_type PetscVector<T>::local_size () const 01088 { 01089 libmesh_assert (this->initialized()); 01090 01091 PetscErrorCode ierr=0; 01092 PetscInt petsc_size=0; 01093 01094 ierr = VecGetLocalSize(_vec, &petsc_size); 01095 LIBMESH_CHKERRABORT(ierr); 01096 01097 return static_cast<numeric_index_type>(petsc_size); 01098 } 01099 01100 01101 01102 template <typename T> 01103 inline 01104 numeric_index_type PetscVector<T>::first_local_index () const 01105 { 01106 libmesh_assert (this->initialized()); 01107 01108 numeric_index_type first = 0; 01109 01110 if(_array_is_present) // Can we use cached values? 01111 first = _first; 01112 else 01113 { 01114 PetscErrorCode ierr=0; 01115 PetscInt petsc_first=0, petsc_last=0; 01116 ierr = VecGetOwnershipRange (_vec, &petsc_first, &petsc_last); 01117 LIBMESH_CHKERRABORT(ierr); 01118 first = static_cast<numeric_index_type>(petsc_first); 01119 } 01120 01121 return first; 01122 } 01123 01124 01125 01126 template <typename T> 01127 inline 01128 numeric_index_type PetscVector<T>::last_local_index () const 01129 { 01130 libmesh_assert (this->initialized()); 01131 01132 numeric_index_type last = 0; 01133 01134 if(_array_is_present) // Can we use cached values? 01135 last = _last; 01136 else 01137 { 01138 PetscErrorCode ierr=0; 01139 PetscInt petsc_first=0, petsc_last=0; 01140 ierr = VecGetOwnershipRange (_vec, &petsc_first, &petsc_last); 01141 LIBMESH_CHKERRABORT(ierr); 01142 last = static_cast<numeric_index_type>(petsc_last); 01143 } 01144 01145 return last; 01146 } 01147 01148 01149 01150 template <typename T> 01151 inline 01152 numeric_index_type PetscVector<T>::map_global_to_local_index (const numeric_index_type i) const 01153 { 01154 libmesh_assert (this->initialized()); 01155 01156 numeric_index_type first=0; 01157 numeric_index_type last=0; 01158 01159 if(_array_is_present) // Can we use cached values? 01160 { 01161 first = _first; 01162 last = _last; 01163 } 01164 else 01165 { 01166 PetscErrorCode ierr=0; 01167 PetscInt petsc_first=0, petsc_last=0; 01168 ierr = VecGetOwnershipRange (_vec, &petsc_first, &petsc_last); 01169 LIBMESH_CHKERRABORT(ierr); 01170 first = static_cast<numeric_index_type>(petsc_first); 01171 last = static_cast<numeric_index_type>(petsc_last); 01172 } 01173 01174 01175 if((i>=first) && (i<last)) 01176 { 01177 return i-first; 01178 } 01179 01180 GlobalToLocalMap::const_iterator it = _global_to_local_map.find(i); 01181 #ifndef NDEBUG 01182 const GlobalToLocalMap::const_iterator end = _global_to_local_map.end(); 01183 if (it == end) 01184 { 01185 std::ostringstream error_message; 01186 error_message << "No index " << i << " in ghosted vector.\n" 01187 << "Vector contains [" << first << ',' << last << ")\n"; 01188 GlobalToLocalMap::const_iterator b = _global_to_local_map.begin(); 01189 if (b == end) 01190 { 01191 error_message << "And empty ghost array.\n"; 01192 } 01193 else 01194 { 01195 error_message << "And ghost array {" << b->first; 01196 for (++b; b != end; ++b) 01197 error_message << ',' << b->first; 01198 error_message << "}\n"; 01199 } 01200 01201 libmesh_error_msg(error_message.str()); 01202 } 01203 libmesh_assert (it != _global_to_local_map.end()); 01204 #endif 01205 return it->second+last-first; 01206 } 01207 01208 01209 01210 template <typename T> 01211 inline 01212 T PetscVector<T>::operator() (const numeric_index_type i) const 01213 { 01214 this->_get_array(); 01215 01216 const numeric_index_type local_index = this->map_global_to_local_index(i); 01217 01218 #ifndef NDEBUG 01219 if(this->type() == GHOSTED) 01220 { 01221 libmesh_assert_less (local_index, _local_size); 01222 } 01223 #endif 01224 01225 return static_cast<T>(_values[local_index]); 01226 } 01227 01228 01229 01230 template <typename T> 01231 inline 01232 void PetscVector<T>::get(const std::vector<numeric_index_type>& index, T* values) const 01233 { 01234 this->_get_array(); 01235 01236 const std::size_t num = index.size(); 01237 01238 for(std::size_t i=0; i<num; i++) 01239 { 01240 const numeric_index_type local_index = this->map_global_to_local_index(index[i]); 01241 #ifndef NDEBUG 01242 if(this->type() == GHOSTED) 01243 { 01244 libmesh_assert_less (local_index, _local_size); 01245 } 01246 #endif 01247 values[i] = static_cast<T>(_values[local_index]); 01248 } 01249 } 01250 01251 01252 01253 template <typename T> 01254 inline 01255 Real PetscVector<T>::min () const 01256 { 01257 this->_restore_array(); 01258 01259 PetscErrorCode ierr=0; 01260 PetscInt index=0; 01261 PetscReal returnval=0.; 01262 01263 ierr = VecMin (_vec, &index, &returnval); 01264 LIBMESH_CHKERRABORT(ierr); 01265 01266 // this return value is correct: VecMin returns a PetscReal 01267 return static_cast<Real>(returnval); 01268 } 01269 01270 01271 01272 template <typename T> 01273 inline 01274 Real PetscVector<T>::max() const 01275 { 01276 this->_restore_array(); 01277 01278 PetscErrorCode ierr=0; 01279 PetscInt index=0; 01280 PetscReal returnval=0.; 01281 01282 ierr = VecMax (_vec, &index, &returnval); 01283 LIBMESH_CHKERRABORT(ierr); 01284 01285 // this return value is correct: VecMax returns a PetscReal 01286 return static_cast<Real>(returnval); 01287 } 01288 01289 01290 01291 template <typename T> 01292 inline 01293 void PetscVector<T>::swap (NumericVector<T> &other) 01294 { 01295 NumericVector<T>::swap(other); 01296 01297 PetscVector<T>& v = cast_ref<PetscVector<T>&>(other); 01298 01299 std::swap(_vec, v._vec); 01300 std::swap(_destroy_vec_on_exit, v._destroy_vec_on_exit); 01301 std::swap(_global_to_local_map, v._global_to_local_map); 01302 std::swap(_array_is_present, v._array_is_present); 01303 std::swap(_local_form, v._local_form); 01304 std::swap(_values, v._values); 01305 } 01306 01307 01308 01309 template <typename T> 01310 inline 01311 void PetscVector<T>::_get_array(void) const 01312 { 01313 libmesh_assert (this->initialized()); 01314 if(!_array_is_present) 01315 { 01316 PetscErrorCode ierr=0; 01317 if(this->type() != GHOSTED) 01318 { 01319 ierr = VecGetArray(_vec, &_values); 01320 LIBMESH_CHKERRABORT(ierr); 01321 } 01322 else 01323 { 01324 ierr = VecGhostGetLocalForm (_vec,&_local_form); 01325 LIBMESH_CHKERRABORT(ierr); 01326 ierr = VecGetArray(_local_form, &_values); 01327 LIBMESH_CHKERRABORT(ierr); 01328 #ifndef NDEBUG 01329 PetscInt my_local_size = 0; 01330 ierr = VecGetLocalSize(_local_form, &my_local_size); 01331 LIBMESH_CHKERRABORT(ierr); 01332 _local_size = static_cast<numeric_index_type>(my_local_size); 01333 #endif 01334 } 01335 01336 { // cache ownership range 01337 PetscInt petsc_first=0, petsc_last=0; 01338 ierr = VecGetOwnershipRange (_vec, &petsc_first, &petsc_last); 01339 LIBMESH_CHKERRABORT(ierr); 01340 _first = static_cast<numeric_index_type>(petsc_first); 01341 _last = static_cast<numeric_index_type>(petsc_last); 01342 } 01343 01344 _array_is_present = true; 01345 } 01346 } 01347 01348 01349 01350 template <typename T> 01351 inline 01352 void PetscVector<T>::_restore_array(void) const 01353 { 01354 libmesh_assert (this->initialized()); 01355 if(_array_is_present) 01356 { 01357 PetscErrorCode ierr=0; 01358 if(this->type() != GHOSTED) 01359 { 01360 ierr = VecRestoreArray (_vec, &_values); 01361 LIBMESH_CHKERRABORT(ierr); 01362 _values = NULL; 01363 } 01364 else 01365 { 01366 ierr = VecRestoreArray (_local_form, &_values); 01367 LIBMESH_CHKERRABORT(ierr); 01368 _values = NULL; 01369 ierr = VecGhostRestoreLocalForm (_vec,&_local_form); 01370 LIBMESH_CHKERRABORT(ierr); 01371 _local_form = NULL; 01372 #ifndef NDEBUG 01373 _local_size = 0; 01374 #endif 01375 } 01376 _array_is_present = false; 01377 } 01378 } 01379 01380 01381 #ifdef LIBMESH_HAVE_CXX11 01382 static_assert(sizeof(PetscInt) == sizeof(numeric_index_type), 01383 "PETSc and libMesh integer sizes must match!"); 01384 #endif 01385 01386 01387 inline 01388 PetscInt* numeric_petsc_cast(const numeric_index_type* p) 01389 { 01390 return reinterpret_cast<PetscInt*>(const_cast<numeric_index_type*>(p)); 01391 } 01392 01393 } // namespace libMesh 01394 01395 01396 #endif // #ifdef LIBMESH_HAVE_PETSC 01397 #endif // LIBMESH_PETSC_VECTOR_H