$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 #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 cast_int<unsigned int>(_val.size()); 00083 } 00084 00088 virtual bool empty() const { return _val.empty(); } 00089 00093 virtual void zero(); 00094 00098 const T & operator() (const unsigned int i) const; 00099 00103 T & operator() (const unsigned int i); 00104 00108 virtual T el(const unsigned int i) const { return (*this)(i); } 00109 00113 virtual T & el(const unsigned int i) { return (*this)(i); } 00114 00118 template <typename T2> 00119 DenseVector<T>& operator = (const DenseVector<T2>& other_vector); 00120 00124 void swap(DenseVector<T>& other_vector); 00125 00129 void resize (const unsigned int n); 00130 00134 void scale (const T factor); 00135 00139 DenseVector<T>& operator*= (const T factor); 00140 00146 template <typename T2, typename T3> 00147 typename boostcopy::enable_if_c< 00148 ScalarTraits<T2>::value, void >::type 00149 add (const T2 factor, 00150 const DenseVector<T3>& vec); 00151 00156 template <typename T2> 00157 typename CompareTypes<T, T2>::supertype dot (const DenseVector<T2> &vec) const; 00158 00163 template <typename T2> 00164 typename CompareTypes<T, T2>::supertype indefinite_dot (const DenseVector<T2> &vec) const; 00165 00169 template <typename T2> 00170 bool operator== (const DenseVector<T2> &vec) const; 00171 00175 template <typename T2> 00176 bool operator!= (const DenseVector<T2> &vec) const; 00177 00181 template <typename T2> 00182 DenseVector<T>& operator+= (const DenseVector<T2> &vec); 00183 00187 template <typename T2> 00188 DenseVector<T>& operator-= (const DenseVector<T2> &vec); 00189 00195 Real min () const; 00196 00202 Real max () const; 00203 00208 Real l1_norm () const; 00209 00215 Real l2_norm () const; 00216 00222 Real linfty_norm () const; 00223 00228 void get_principal_subvector (unsigned int sub_n, DenseVector<T>& dest) const; 00229 00235 std::vector<T>& get_values() { return _val; } 00236 00242 const std::vector<T>& get_values() const { return _val; } 00243 00244 private: 00245 00249 std::vector<T> _val; 00250 00251 }; 00252 00253 00254 00255 // ------------------------------------------------------------ 00256 // DenseVector member functions 00257 template<typename T> 00258 inline 00259 DenseVector<T>::DenseVector(const unsigned int n) : 00260 _val (n, T(0.)) 00261 { 00262 } 00263 00264 00265 00266 template<typename T> 00267 template<typename T2> 00268 inline 00269 DenseVector<T>::DenseVector (const DenseVector<T2>& other_vector) : 00270 DenseVectorBase<T>() 00271 { 00272 const std::vector<T2> &other_vals = other_vector.get_values(); 00273 00274 _val.clear(); 00275 _val.reserve(other_vals.size()); 00276 00277 for (unsigned int i=0; i<other_vals.size(); i++) 00278 _val.push_back(other_vals[i]); 00279 } 00280 00281 00282 00283 template<typename T> 00284 template<typename T2> 00285 inline 00286 DenseVector<T>::DenseVector (const std::vector<T2>& other_vector) : 00287 _val(other_vector) 00288 { 00289 } 00290 00291 00292 00293 00294 00295 template<typename T> 00296 template<typename T2> 00297 inline 00298 DenseVector<T>& DenseVector<T>::operator = (const DenseVector<T2>& other_vector) 00299 { 00300 // _val = other_vector._val; 00301 00302 const std::vector<T2> &other_vals = other_vector.get_values(); 00303 00304 _val.clear(); 00305 _val.reserve(other_vals.size()); 00306 00307 for (unsigned int i=0; i<other_vals.size(); i++) 00308 _val.push_back(other_vals[i]); 00309 00310 return *this; 00311 } 00312 00313 00314 00315 template<typename T> 00316 inline 00317 void DenseVector<T>::swap(DenseVector<T>& other_vector) 00318 { 00319 _val.swap(other_vector._val); 00320 } 00321 00322 00323 00324 template<typename T> 00325 inline 00326 void DenseVector<T>::resize(const unsigned int n) 00327 { 00328 _val.resize(n); 00329 00330 zero(); 00331 } 00332 00333 00334 00335 template<typename T> 00336 inline 00337 void DenseVector<T>::zero() 00338 { 00339 std::fill (_val.begin(), 00340 _val.end(), 00341 T(0.)); 00342 } 00343 00344 00345 00346 template<typename T> 00347 inline 00348 const T & DenseVector<T>::operator () (const unsigned int i) const 00349 { 00350 libmesh_assert_less (i, _val.size()); 00351 00352 return _val[i]; 00353 } 00354 00355 00356 00357 template<typename T> 00358 inline 00359 T & DenseVector<T>::operator () (const unsigned int i) 00360 { 00361 libmesh_assert_less (i, _val.size()); 00362 00363 return _val[i]; 00364 } 00365 00366 00367 00368 template<typename T> 00369 inline 00370 void DenseVector<T>::scale (const T factor) 00371 { 00372 for (unsigned int i=0; i<_val.size(); i++) 00373 _val[i] *= factor; 00374 } 00375 00376 00377 00378 template<typename T> 00379 inline 00380 DenseVector<T>& DenseVector<T>::operator*= (const T factor) 00381 { 00382 this->scale(factor); 00383 return *this; 00384 } 00385 00386 00387 00388 template<typename T> 00389 template<typename T2, typename T3> 00390 inline 00391 typename boostcopy::enable_if_c< 00392 ScalarTraits<T2>::value, void >::type 00393 DenseVector<T>::add (const T2 factor, 00394 const DenseVector<T3>& vec) 00395 { 00396 libmesh_assert_equal_to (this->size(), vec.size()); 00397 00398 for (unsigned int i=0; i<this->size(); i++) 00399 (*this)(i) += static_cast<T>(factor)*vec(i); 00400 } 00401 00402 template<typename T> 00403 template<typename T2> 00404 inline 00405 typename CompareTypes<T, T2>::supertype DenseVector<T>::dot (const DenseVector<T2>& vec) const 00406 { 00407 libmesh_assert_equal_to (this->size(), vec.size()); 00408 00409 typename CompareTypes<T, T2>::supertype val = 0.; 00410 00411 for (unsigned int i=0; i<this->size(); i++) 00412 val += (*this)(i)*libmesh_conj(vec(i)); 00413 00414 return val; 00415 } 00416 00417 template<typename T> 00418 template<typename T2> 00419 inline 00420 typename CompareTypes<T, T2>::supertype DenseVector<T>::indefinite_dot (const DenseVector<T2>& vec) const 00421 { 00422 libmesh_assert_equal_to (this->size(), vec.size()); 00423 00424 typename CompareTypes<T, T2>::supertype val = 0.; 00425 00426 for (unsigned int i=0; i<this->size(); i++) 00427 val += (*this)(i)*(vec(i)); 00428 00429 return val; 00430 } 00431 00432 template<typename T> 00433 template<typename T2> 00434 inline 00435 bool DenseVector<T>::operator== (const DenseVector<T2>& vec) const 00436 { 00437 libmesh_assert_equal_to (this->size(), vec.size()); 00438 00439 for (unsigned int i=0; i<this->size(); i++) 00440 if ((*this)(i) != vec(i)) 00441 return false; 00442 00443 return true; 00444 } 00445 00446 00447 00448 template<typename T> 00449 template<typename T2> 00450 inline 00451 bool DenseVector<T>::operator!= (const DenseVector<T2>& vec) const 00452 { 00453 libmesh_assert_equal_to (this->size(), vec.size()); 00454 00455 for (unsigned int i=0; i<this->size(); i++) 00456 if ((*this)(i) != vec(i)) 00457 return true; 00458 00459 return false; 00460 } 00461 00462 00463 00464 template<typename T> 00465 template<typename T2> 00466 inline 00467 DenseVector<T>& DenseVector<T>::operator+= (const DenseVector<T2>& vec) 00468 { 00469 libmesh_assert_equal_to (this->size(), vec.size()); 00470 00471 for (unsigned int i=0; i<this->size(); i++) 00472 (*this)(i) += vec(i); 00473 00474 return *this; 00475 } 00476 00477 00478 00479 template<typename T> 00480 template<typename T2> 00481 inline 00482 DenseVector<T>& DenseVector<T>::operator-= (const DenseVector<T2>& vec) 00483 { 00484 libmesh_assert_equal_to (this->size(), vec.size()); 00485 00486 for (unsigned int i=0; i<this->size(); i++) 00487 (*this)(i) -= vec(i); 00488 00489 return *this; 00490 } 00491 00492 00493 00494 template<typename T> 00495 inline 00496 Real DenseVector<T>::min () const 00497 { 00498 libmesh_assert (this->size()); 00499 Real my_min = libmesh_real((*this)(0)); 00500 00501 for (unsigned int i=1; i!=this->size(); i++) 00502 { 00503 Real current = libmesh_real((*this)(i)); 00504 my_min = (my_min < current? my_min : current); 00505 } 00506 return my_min; 00507 } 00508 00509 00510 00511 template<typename T> 00512 inline 00513 Real DenseVector<T>::max () const 00514 { 00515 libmesh_assert (this->size()); 00516 Real my_max = libmesh_real((*this)(0)); 00517 00518 for (unsigned int i=1; i!=this->size(); i++) 00519 { 00520 Real current = libmesh_real((*this)(i)); 00521 my_max = (my_max > current? my_max : current); 00522 } 00523 return my_max; 00524 } 00525 00526 00527 00528 template<typename T> 00529 inline 00530 Real DenseVector<T>::l1_norm () const 00531 { 00532 Real my_norm = 0.; 00533 for (unsigned int i=0; i!=this->size(); i++) 00534 { 00535 my_norm += std::abs((*this)(i)); 00536 } 00537 return my_norm; 00538 } 00539 00540 00541 00542 template<typename T> 00543 inline 00544 Real DenseVector<T>::l2_norm () const 00545 { 00546 Real my_norm = 0.; 00547 for (unsigned int i=0; i!=this->size(); i++) 00548 { 00549 my_norm += TensorTools::norm_sq((*this)(i)); 00550 } 00551 return sqrt(my_norm); 00552 } 00553 00554 00555 00556 template<typename T> 00557 inline 00558 Real DenseVector<T>::linfty_norm () const 00559 { 00560 if (!this->size()) 00561 return 0.; 00562 Real my_norm = TensorTools::norm_sq((*this)(0)); 00563 00564 for (unsigned int i=1; i!=this->size(); i++) 00565 { 00566 Real current = TensorTools::norm_sq((*this)(i)); 00567 my_norm = (my_norm > current? my_norm : current); 00568 } 00569 return sqrt(my_norm); 00570 } 00571 00572 template<typename T> 00573 inline 00574 void DenseVector<T>::get_principal_subvector (unsigned int sub_n, 00575 DenseVector<T>& dest) const 00576 { 00577 libmesh_assert_less_equal ( sub_n, this->size() ); 00578 00579 dest.resize(sub_n); 00580 for(unsigned int i=0; i<sub_n; i++) 00581 dest(i) = (*this)(i); 00582 } 00583 00584 } // namespace libMesh 00585 00586 #endif // LIBMESH_DENSE_VECTOR_H