dense_vector.h
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 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 #ifndef LIBMESH_DENSE_VECTOR_H 00021 #define LIBMESH_DENSE_VECTOR_H 00022 00023 // Local Includes 00024 #include "libmesh/libmesh_common.h" 00025 #include "libmesh/compare_types.h" 00026 #include "libmesh/dense_vector_base.h" 00027 #include "libmesh/tensor_tools.h" 00028 00029 // C++ includes 00030 #include <vector> 00031 00032 namespace libMesh 00033 { 00034 00035 // Forward Declarations 00036 00037 00038 00048 // ------------------------------------------------------------ 00049 // DenseVector class definition 00050 template<typename T> 00051 class DenseVector : public DenseVectorBase<T> 00052 { 00053 public: 00054 00058 explicit 00059 DenseVector(const unsigned int n=0); 00060 00064 template <typename T2> 00065 DenseVector (const DenseVector<T2>& other_vector); 00066 00070 template <typename T2> 00071 DenseVector (const std::vector<T2>& other_vector); 00072 00076 ~DenseVector() {} 00077 00081 virtual unsigned int size() const { 00082 return libmesh_cast_int<unsigned int>(_val.size()); 00083 } 00084 00088 virtual void zero(); 00089 00093 T operator() (const unsigned int i) const; 00094 00098 T & operator() (const unsigned int i); 00099 00103 virtual T el(const unsigned int i) const { return (*this)(i); } 00104 00108 virtual T & el(const unsigned int i) { return (*this)(i); } 00109 00113 template <typename T2> 00114 DenseVector<T>& operator = (const DenseVector<T2>& other_vector); 00115 00119 void swap(DenseVector<T>& other_vector); 00120 00124 void resize (const unsigned int n); 00125 00129 void scale (const T factor); 00130 00134 DenseVector<T>& operator*= (const T factor); 00135 00141 template <typename T2, typename T3> 00142 typename boostcopy::enable_if_c< 00143 ScalarTraits<T2>::value, void >::type 00144 add (const T2 factor, 00145 const DenseVector<T3>& vec); 00146 00151 template <typename T2> 00152 Number dot (const DenseVector<T2> &vec) const; 00153 00158 template <typename T2> 00159 Number indefinite_dot (const DenseVector<T2> &vec) const; 00160 00164 template <typename T2> 00165 bool operator== (const DenseVector<T2> &vec) const; 00166 00170 template <typename T2> 00171 bool operator!= (const DenseVector<T2> &vec) const; 00172 00176 template <typename T2> 00177 DenseVector<T>& operator+= (const DenseVector<T2> &vec); 00178 00182 template <typename T2> 00183 DenseVector<T>& operator-= (const DenseVector<T2> &vec); 00184 00190 Real min () const; 00191 00197 Real max () const; 00198 00203 Real l1_norm () const; 00204 00210 Real l2_norm () const; 00211 00217 Real linfty_norm () const; 00218 00223 void get_principal_subvector (unsigned int sub_n, DenseVector<T>& dest) const; 00224 00230 std::vector<T>& get_values() { return _val; } 00231 00237 const std::vector<T>& get_values() const { return _val; } 00238 00239 private: 00240 00244 std::vector<T> _val; 00245 00246 }; 00247 00248 00249 00250 // ------------------------------------------------------------ 00251 // DenseVector member functions 00252 template<typename T> 00253 inline 00254 DenseVector<T>::DenseVector(const unsigned int n) : 00255 _val (n, T(0.)) 00256 { 00257 } 00258 00259 00260 00261 template<typename T> 00262 template<typename T2> 00263 inline 00264 DenseVector<T>::DenseVector (const DenseVector<T2>& other_vector) : 00265 DenseVectorBase<T>() 00266 { 00267 const std::vector<T2> &other_vals = other_vector.get_values(); 00268 00269 _val.clear(); 00270 _val.reserve(other_vals.size()); 00271 00272 for (unsigned int i=0; i<other_vals.size(); i++) 00273 _val.push_back(other_vals[i]); 00274 } 00275 00276 00277 00278 template<typename T> 00279 template<typename T2> 00280 inline 00281 DenseVector<T>::DenseVector (const std::vector<T2>& other_vector) : 00282 _val(other_vector) 00283 { 00284 } 00285 00286 00287 00288 00289 00290 template<typename T> 00291 template<typename T2> 00292 inline 00293 DenseVector<T>& DenseVector<T>::operator = (const DenseVector<T2>& other_vector) 00294 { 00295 // _val = other_vector._val; 00296 00297 const std::vector<T2> &other_vals = other_vector.get_values(); 00298 00299 _val.clear(); 00300 _val.reserve(other_vals.size()); 00301 00302 for (unsigned int i=0; i<other_vals.size(); i++) 00303 _val.push_back(other_vals[i]); 00304 00305 return *this; 00306 } 00307 00308 00309 00310 template<typename T> 00311 inline 00312 void DenseVector<T>::swap(DenseVector<T>& other_vector) 00313 { 00314 _val.swap(other_vector._val); 00315 } 00316 00317 00318 00319 template<typename T> 00320 inline 00321 void DenseVector<T>::resize(const unsigned int n) 00322 { 00323 _val.resize(n); 00324 00325 zero(); 00326 } 00327 00328 00329 00330 template<typename T> 00331 inline 00332 void DenseVector<T>::zero() 00333 { 00334 std::fill (_val.begin(), 00335 _val.end(), 00336 T(0.)); 00337 } 00338 00339 00340 00341 template<typename T> 00342 inline 00343 T DenseVector<T>::operator () (const unsigned int i) const 00344 { 00345 libmesh_assert_less (i, _val.size()); 00346 00347 return _val[i]; 00348 } 00349 00350 00351 00352 template<typename T> 00353 inline 00354 T & DenseVector<T>::operator () (const unsigned int i) 00355 { 00356 libmesh_assert_less (i, _val.size()); 00357 00358 return _val[i]; 00359 } 00360 00361 00362 00363 template<typename T> 00364 inline 00365 void DenseVector<T>::scale (const T factor) 00366 { 00367 for (unsigned int i=0; i<_val.size(); i++) 00368 _val[i] *= factor; 00369 } 00370 00371 00372 00373 template<typename T> 00374 inline 00375 DenseVector<T>& DenseVector<T>::operator*= (const T factor) 00376 { 00377 this->scale(factor); 00378 return *this; 00379 } 00380 00381 00382 00383 template<typename T> 00384 template<typename T2, typename T3> 00385 inline 00386 typename boostcopy::enable_if_c< 00387 ScalarTraits<T2>::value, void >::type 00388 DenseVector<T>::add (const T2 factor, 00389 const DenseVector<T3>& vec) 00390 { 00391 libmesh_assert_equal_to (this->size(), vec.size()); 00392 00393 for (unsigned int i=0; i<this->size(); i++) 00394 (*this)(i) += factor*vec(i); 00395 } 00396 00397 template<typename T> 00398 template<typename T2> 00399 inline 00400 Number DenseVector<T>::dot (const DenseVector<T2>& vec) const 00401 { 00402 libmesh_assert_equal_to (this->size(), vec.size()); 00403 00404 Number val = 0.; 00405 00406 for (unsigned int i=0; i<this->size(); i++) 00407 val += (*this)(i)*libmesh_conj(vec(i)); 00408 00409 return val; 00410 } 00411 00412 template<typename T> 00413 template<typename T2> 00414 inline 00415 Number DenseVector<T>::indefinite_dot (const DenseVector<T2>& vec) const 00416 { 00417 libmesh_assert_equal_to (this->size(), vec.size()); 00418 00419 Number val = 0.; 00420 00421 for (unsigned int i=0; i<this->size(); i++) 00422 val += (*this)(i)*(vec(i)); 00423 00424 return val; 00425 } 00426 00427 template<typename T> 00428 template<typename T2> 00429 inline 00430 bool DenseVector<T>::operator== (const DenseVector<T2>& vec) const 00431 { 00432 libmesh_assert_equal_to (this->size(), vec.size()); 00433 00434 for (unsigned int i=0; i<this->size(); i++) 00435 if ((*this)(i) != vec(i)) 00436 return false; 00437 00438 return true; 00439 } 00440 00441 00442 00443 template<typename T> 00444 template<typename T2> 00445 inline 00446 bool DenseVector<T>::operator!= (const DenseVector<T2>& vec) const 00447 { 00448 libmesh_assert_equal_to (this->size(), vec.size()); 00449 00450 for (unsigned int i=0; i<this->size(); i++) 00451 if ((*this)(i) != vec(i)) 00452 return true; 00453 00454 return false; 00455 } 00456 00457 00458 00459 template<typename T> 00460 template<typename T2> 00461 inline 00462 DenseVector<T>& DenseVector<T>::operator+= (const DenseVector<T2>& vec) 00463 { 00464 libmesh_assert_equal_to (this->size(), vec.size()); 00465 00466 for (unsigned int i=0; i<this->size(); i++) 00467 (*this)(i) += vec(i); 00468 00469 return *this; 00470 } 00471 00472 00473 00474 template<typename T> 00475 template<typename T2> 00476 inline 00477 DenseVector<T>& DenseVector<T>::operator-= (const DenseVector<T2>& vec) 00478 { 00479 libmesh_assert_equal_to (this->size(), vec.size()); 00480 00481 for (unsigned int i=0; i<this->size(); i++) 00482 (*this)(i) -= vec(i); 00483 00484 return *this; 00485 } 00486 00487 00488 00489 template<typename T> 00490 inline 00491 Real DenseVector<T>::min () const 00492 { 00493 libmesh_assert (this->size()); 00494 Real my_min = libmesh_real((*this)(0)); 00495 00496 for (unsigned int i=1; i!=this->size(); i++) 00497 { 00498 Real current = libmesh_real((*this)(i)); 00499 my_min = (my_min < current? my_min : current); 00500 } 00501 return my_min; 00502 } 00503 00504 00505 00506 template<typename T> 00507 inline 00508 Real DenseVector<T>::max () const 00509 { 00510 libmesh_assert (this->size()); 00511 Real my_max = libmesh_real((*this)(0)); 00512 00513 for (unsigned int i=1; i!=this->size(); i++) 00514 { 00515 Real current = libmesh_real((*this)(i)); 00516 my_max = (my_max > current? my_max : current); 00517 } 00518 return my_max; 00519 } 00520 00521 00522 00523 template<typename T> 00524 inline 00525 Real DenseVector<T>::l1_norm () const 00526 { 00527 Real my_norm = 0.; 00528 for (unsigned int i=0; i!=this->size(); i++) 00529 { 00530 my_norm += std::abs((*this)(i)); 00531 } 00532 return my_norm; 00533 } 00534 00535 00536 00537 template<typename T> 00538 inline 00539 Real DenseVector<T>::l2_norm () const 00540 { 00541 Real my_norm = 0.; 00542 for (unsigned int i=0; i!=this->size(); i++) 00543 { 00544 my_norm += TensorTools::norm_sq((*this)(i)); 00545 } 00546 return sqrt(my_norm); 00547 } 00548 00549 00550 00551 template<typename T> 00552 inline 00553 Real DenseVector<T>::linfty_norm () const 00554 { 00555 if (!this->size()) 00556 return 0.; 00557 Real my_norm = TensorTools::norm_sq((*this)(0)); 00558 00559 for (unsigned int i=1; i!=this->size(); i++) 00560 { 00561 Real current = TensorTools::norm_sq((*this)(i)); 00562 my_norm = (my_norm > current? my_norm : current); 00563 } 00564 return sqrt(my_norm); 00565 } 00566 00567 template<typename T> 00568 inline 00569 void DenseVector<T>::get_principal_subvector (unsigned int sub_n, 00570 DenseVector<T>& dest) const 00571 { 00572 libmesh_assert_less_equal ( sub_n, this->size() ); 00573 00574 dest.resize(sub_n); 00575 for(unsigned int i=0; i<sub_n; i++) 00576 dest(i) = (*this)(i); 00577 } 00578 00579 } // namespace libMesh 00580 00581 #endif // LIBMESH_DENSE_VECTOR_H 00582
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: