$extrastylesheet
eigen_sparse_vector.C
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 // C++ includes
00021 #include <algorithm> // for std::min
00022 #include <limits>
00023 
00024 // Local Includes
00025 #include "libmesh/dense_subvector.h"
00026 #include "libmesh/dense_vector.h"
00027 #include "libmesh/eigen_sparse_vector.h"
00028 #include "libmesh/eigen_sparse_matrix.h"
00029 
00030 
00031 #ifdef LIBMESH_HAVE_EIGEN
00032 
00033 namespace libMesh
00034 {
00035 
00036 template <typename T>
00037 T EigenSparseVector<T>::sum () const
00038 {
00039   libmesh_assert (this->closed());
00040   libmesh_assert (this->initialized());
00041 
00042   return _vec.sum();
00043 }
00044 
00045 
00046 
00047 template <typename T>
00048 Real EigenSparseVector<T>::l1_norm () const
00049 {
00050   libmesh_assert (this->closed());
00051   libmesh_assert (this->initialized());
00052 
00053   return _vec.lpNorm<1>();
00054 }
00055 
00056 
00057 
00058 template <typename T>
00059 Real EigenSparseVector<T>::l2_norm () const
00060 {
00061   libmesh_assert (this->closed());
00062   libmesh_assert (this->initialized());
00063 
00064   return _vec.lpNorm<2>();
00065 }
00066 
00067 
00068 
00069 template <typename T>
00070 Real EigenSparseVector<T>::linfty_norm () const
00071 {
00072   libmesh_assert (this->closed());
00073   libmesh_assert (this->initialized());
00074 
00075   return _vec.lpNorm<Eigen::Infinity>();
00076 }
00077 
00078 
00079 
00080 template <typename T>
00081 NumericVector<T>& EigenSparseVector<T>::operator += (const NumericVector<T>& v_in)
00082 {
00083   libmesh_assert (this->closed());
00084 
00085   const EigenSparseVector<T>& v = cast_ref<const EigenSparseVector<T>&>(v_in);
00086 
00087   _vec += v._vec;
00088 
00089   return *this;
00090 }
00091 
00092 
00093 
00094 
00095 template <typename T>
00096 NumericVector<T>& EigenSparseVector<T>::operator -= (const NumericVector<T>& v_in)
00097 {
00098   libmesh_assert (this->closed());
00099 
00100   const EigenSparseVector<T>& v = cast_ref<const EigenSparseVector<T>&>(v_in);
00101 
00102   _vec -= v._vec;
00103 
00104   return *this;
00105 }
00106 
00107 
00108 
00109 template <typename T>
00110 NumericVector<T> & EigenSparseVector<T>::operator /= (NumericVector<T> & v_in)
00111 {
00112   libmesh_assert (this->closed());
00113   libmesh_assert_equal_to(size(), v_in.size());
00114 
00115   const EigenSparseVector<T>& v = cast_ref<const EigenSparseVector<T>&>(v_in);
00116 
00117   _vec = _vec.cwiseQuotient(v._vec);
00118 
00119   return *this;
00120 }
00121 
00122 
00123 
00124 
00125 template <typename T>
00126 void EigenSparseVector<T>::reciprocal()
00127 {
00128 #ifndef NDEBUG
00129   const numeric_index_type n = this->size();
00130 
00131   for (numeric_index_type i=0; i<n; i++)
00132     // Don't divide by zero!
00133     libmesh_assert_not_equal_to ((*this)(i), T(0));
00134 #endif
00135 
00136   _vec = _vec.cwiseInverse();
00137 }
00138 
00139 
00140 
00141 template <typename T>
00142 void EigenSparseVector<T>::conjugate()
00143 {
00144   _vec = _vec.conjugate();
00145 }
00146 
00147 
00148 
00149 template <typename T>
00150 void EigenSparseVector<T>::add (const T v)
00151 {
00152   _vec += EigenSV::Constant(this->size(), v);
00153 
00154 #ifndef NDEBUG
00155   this->_is_closed = false;
00156 #endif
00157 }
00158 
00159 
00160 
00161 
00162 template <typename T>
00163 void EigenSparseVector<T>::add (const NumericVector<T>& v_in)
00164 {
00165   libmesh_assert (this->initialized());
00166 
00167   const EigenSparseVector<T>& v = cast_ref<const EigenSparseVector<T>&>(v_in);
00168 
00169   _vec += v._vec;
00170 }
00171 
00172 
00173 
00174 template <typename T>
00175 void EigenSparseVector<T>::add (const T a, const NumericVector<T>& v_in)
00176 {
00177   libmesh_assert (this->initialized());
00178 
00179   const EigenSparseVector<T>& v = cast_ref<const EigenSparseVector<T>&>(v_in);
00180 
00181   _vec += v._vec*a;
00182 }
00183 
00184 
00185 
00186 template <typename T>
00187 void EigenSparseVector<T>::add_vector (const NumericVector<T> &vec_in,
00188                                        const SparseMatrix<T>  &mat_in)
00189 {
00190   // Make sure the data passed in are really in Eigen types
00191   const EigenSparseVector<T>* e_vec = cast_ptr<const EigenSparseVector<T>*>(&vec_in);
00192   const EigenSparseMatrix<T>* mat = cast_ptr<const EigenSparseMatrix<T>*>(&mat_in);
00193 
00194   libmesh_assert(e_vec);
00195   libmesh_assert(mat);
00196 
00197   _vec += mat->_mat*e_vec->_vec;
00198 }
00199 
00200 
00201 
00202 template <typename T>
00203 void EigenSparseVector<T>::add_vector_transpose (const NumericVector<T> &vec_in,
00204                                                  const SparseMatrix<T>  &mat_in)
00205 {
00206   // Make sure the data passed in are really in Eigen types
00207   const EigenSparseVector<T>* e_vec = cast_ptr<const EigenSparseVector<T>*>(&vec_in);
00208   const EigenSparseMatrix<T>* mat = cast_ptr<const EigenSparseMatrix<T>*>(&mat_in);
00209 
00210   libmesh_assert(e_vec);
00211   libmesh_assert(mat);
00212 
00213   _vec += mat->_mat.transpose()*e_vec->_vec;
00214 }
00215 
00216 
00217 
00218 template <typename T>
00219 void EigenSparseVector<T>::scale (const T factor)
00220 {
00221   libmesh_assert (this->initialized());
00222 
00223   _vec *= factor;
00224 }
00225 
00226 
00227 
00228 template <typename T>
00229 void EigenSparseVector<T>::abs()
00230 {
00231   libmesh_assert (this->initialized());
00232 
00233   const numeric_index_type n = this->size();
00234 
00235   for (numeric_index_type i=0; i!=n; ++i)
00236     this->set(i,std::abs((*this)(i)));
00237 }
00238 
00239 
00240 
00241 template <typename T>
00242 T EigenSparseVector<T>::dot (const NumericVector<T>& V) const
00243 {
00244   libmesh_assert (this->initialized());
00245 
00246   // Make sure the NumericVector passed in is really a EigenSparseVector
00247   const EigenSparseVector<T>* v = cast_ptr<const EigenSparseVector<T>*>(&V);
00248   libmesh_assert(v);
00249 
00250   return _vec.dot(v->_vec);
00251 }
00252 
00253 
00254 
00255 template <typename T>
00256 NumericVector<T>&
00257 EigenSparseVector<T>::operator = (const T s)
00258 {
00259   libmesh_assert (this->initialized());
00260   libmesh_assert (this->closed());
00261 
00262   _vec.fill(s);
00263 
00264   return *this;
00265 }
00266 
00267 
00268 
00269 template <typename T>
00270 NumericVector<T>&
00271 EigenSparseVector<T>::operator = (const NumericVector<T>& v_in)
00272 {
00273   // Make sure the NumericVector passed in is really a EigenSparseVector
00274   const EigenSparseVector<T>* v =
00275     cast_ptr<const EigenSparseVector<T>*>(&v_in);
00276 
00277   libmesh_assert(v);
00278 
00279   *this = *v;
00280 
00281   return *this;
00282 }
00283 
00284 
00285 
00286 template <typename T>
00287 EigenSparseVector<T>&
00288 EigenSparseVector<T>::operator = (const EigenSparseVector<T>& v)
00289 {
00290   libmesh_assert (this->initialized());
00291   libmesh_assert (v.closed());
00292   libmesh_assert_equal_to (this->size(), v.size());
00293 
00294   _vec = v._vec;
00295 
00296 #ifndef NDEBUG
00297   this->_is_closed = true;
00298 #endif
00299 
00300   return *this;
00301 }
00302 
00303 
00304 
00305 template <typename T>
00306 NumericVector<T>&
00307 EigenSparseVector<T>::operator = (const std::vector<T>& v)
00308 {
00313   if (this->size() == v.size())
00314     for (numeric_index_type i=0; i<v.size(); i++)
00315       this->set (i, v[i]);
00316 
00317   else
00318     libmesh_error_msg("this->size() = " << this->size() << " must be equal to v.size() = " << v.size());
00319 
00320   return *this;
00321 }
00322 
00323 
00324 template <typename T>
00325 void EigenSparseVector<T>::localize (NumericVector<T>& v_local_in) const
00326 {
00327   // Make sure the NumericVector passed in is really a EigenSparseVector
00328   EigenSparseVector<T>* v_local =
00329     cast_ptr<EigenSparseVector<T>*>(&v_local_in);
00330 
00331   libmesh_assert(v_local);
00332 
00333   *v_local = *this;
00334 }
00335 
00336 
00337 
00338 template <typename T>
00339 void EigenSparseVector<T>::localize (NumericVector<T>& v_local_in,
00340                                      const std::vector<numeric_index_type>& libmesh_dbg_var(send_list)) const
00341 {
00342   // Make sure the NumericVector passed in is really a EigenSparseVector
00343   EigenSparseVector<T>* v_local =
00344     cast_ptr<EigenSparseVector<T>*>(&v_local_in);
00345 
00346   libmesh_assert(v_local);
00347   libmesh_assert_less_equal (send_list.size(), v_local->size());
00348 
00349   *v_local = *this;
00350 }
00351 
00352 
00353 
00354 template <typename T>
00355 void EigenSparseVector<T>::localize (const numeric_index_type libmesh_dbg_var(first_local_idx),
00356                                      const numeric_index_type libmesh_dbg_var(last_local_idx),
00357                                      const std::vector<numeric_index_type>& libmesh_dbg_var(send_list))
00358 {
00359   libmesh_assert_equal_to (first_local_idx, 0);
00360   libmesh_assert_equal_to (last_local_idx+1, this->size());
00361 
00362   libmesh_assert_less_equal (send_list.size(), this->size());
00363 
00364 #ifndef NDEBUG
00365   this->_is_closed = true;
00366 #endif
00367 }
00368 
00369 
00370 
00371 template <typename T>
00372 void EigenSparseVector<T>::localize (std::vector<T>& v_local) const
00373 
00374 {
00375   v_local.resize(this->size());
00376 
00377   for (numeric_index_type i=0; i<v_local.size(); i++)
00378     v_local[i] = (*this)(i);
00379 }
00380 
00381 
00382 
00383 template <typename T>
00384 void EigenSparseVector<T>::localize_to_one (std::vector<T>& v_local,
00385                                             const processor_id_type libmesh_dbg_var(pid)) const
00386 {
00387   libmesh_assert_equal_to (pid, 0);
00388 
00389   this->localize (v_local);
00390 }
00391 
00392 
00393 
00394 template <typename T>
00395 void EigenSparseVector<T>::pointwise_mult (const NumericVector<T>& /*vec1*/,
00396                                            const NumericVector<T>& /*vec2*/)
00397 {
00398   libmesh_not_implemented();
00399 }
00400 
00401 
00402 
00403 template <typename T>
00404 Real EigenSparseVector<T>::max() const
00405 {
00406   libmesh_assert (this->initialized());
00407   if (!this->size())
00408     return -std::numeric_limits<Real>::max();
00409 
00410 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
00411   Real the_max = libmesh_real((*this)(0));
00412 
00413   const numeric_index_type n = this->size();
00414 
00415   for (numeric_index_type i=1; i<n; i++)
00416     the_max = std::max (the_max, libmesh_real((*this)(i)));
00417 
00418   return the_max;
00419 #else
00420   return libmesh_real(_vec.maxCoeff());
00421 #endif
00422 }
00423 
00424 
00425 
00426 template <typename T>
00427 Real EigenSparseVector<T>::min () const
00428 {
00429   libmesh_assert (this->initialized());
00430   if (!this->size())
00431     return std::numeric_limits<Real>::max();
00432 
00433 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
00434   Real the_min = libmesh_real((*this)(0));
00435 
00436   const numeric_index_type n = this->size();
00437 
00438   for (numeric_index_type i=1; i<n; i++)
00439     the_min = std::min (the_min, libmesh_real((*this)(i)));
00440 
00441   return the_min;
00442 #else
00443   return libmesh_real(_vec.minCoeff());
00444 #endif
00445 }
00446 
00447 
00448 //------------------------------------------------------------------
00449 // Explicit instantiations
00450 template class EigenSparseVector<Number>;
00451 
00452 } // namespace libMesh
00453 
00454 
00455 #endif // #ifdef LIBMESH_HAVE_EIGEN